xref: /petsc/src/sys/objects/aoptions.c (revision 126553251c28899133a26d03260141881bdb71df)
17d0a6c19SBarry Smith 
21a1499c8SBarry Smith 
353acd3b1SBarry Smith /*
43fc1eb6aSBarry Smith    Implements the higher-level options database querying methods. These are self-documenting and can attach at runtime to
53fc1eb6aSBarry Smith    GUI code to display the options and get values from the users.
653acd3b1SBarry Smith 
753acd3b1SBarry Smith */
853acd3b1SBarry Smith 
9afcb2eb5SJed Brown #include <petsc-private/petscimpl.h>        /*I  "petscsys.h"   I*/
10665c2dedSJed Brown #include <petscviewer.h>
1153acd3b1SBarry Smith 
122aa6d131SJed Brown #define ManSection(str) ((str) ? (str) : "None")
132aa6d131SJed Brown 
1453acd3b1SBarry Smith /*
1553acd3b1SBarry Smith     Keep a linked list of options that have been posted and we are waiting for
163fc1eb6aSBarry Smith    user selection. See the manual page for PetscOptionsBegin()
1753acd3b1SBarry Smith 
1853acd3b1SBarry Smith     Eventually we'll attach this beast to a MPI_Comm
1953acd3b1SBarry Smith */
20e55864a3SBarry Smith 
2153acd3b1SBarry Smith 
2253acd3b1SBarry Smith #undef __FUNCT__
2353acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsBegin_Private"
2453acd3b1SBarry Smith /*
2553acd3b1SBarry Smith     Handles setting up the data structure in a call to PetscOptionsBegin()
2653acd3b1SBarry Smith */
278c34d3f5SBarry Smith PetscErrorCode PetscOptionsBegin_Private(PetscOptions *PetscOptionsObject,MPI_Comm comm,const char prefix[],const char title[],const char mansec[])
2853acd3b1SBarry Smith {
2953acd3b1SBarry Smith   PetscErrorCode ierr;
3053acd3b1SBarry Smith 
3153acd3b1SBarry Smith   PetscFunctionBegin;
32e55864a3SBarry Smith   PetscOptionsObject->next          = 0;
33e55864a3SBarry Smith   PetscOptionsObject->comm          = comm;
34e55864a3SBarry Smith   PetscOptionsObject->changedmethod = PETSC_FALSE;
35a297a907SKarl Rupp 
36e55864a3SBarry Smith   ierr = PetscStrallocpy(prefix,&PetscOptionsObject->prefix);CHKERRQ(ierr);
37e55864a3SBarry Smith   ierr = PetscStrallocpy(title,&PetscOptionsObject->title);CHKERRQ(ierr);
3853acd3b1SBarry Smith 
39e55864a3SBarry Smith   ierr = PetscOptionsHasName(NULL,"-help",&PetscOptionsObject->printhelp);CHKERRQ(ierr);
40e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1) {
41e55864a3SBarry Smith     if (!PetscOptionsObject->alreadyprinted) {
4253acd3b1SBarry Smith       ierr = (*PetscHelpPrintf)(comm,"%s -------------------------------------------------\n",title);CHKERRQ(ierr);
4353acd3b1SBarry Smith     }
4461b37b28SSatish Balay   }
4553acd3b1SBarry Smith   PetscFunctionReturn(0);
4653acd3b1SBarry Smith }
4753acd3b1SBarry Smith 
483194b578SJed Brown #undef __FUNCT__
493194b578SJed Brown #define __FUNCT__ "PetscObjectOptionsBegin_Private"
503194b578SJed Brown /*
513194b578SJed Brown     Handles setting up the data structure in a call to PetscObjectOptionsBegin()
523194b578SJed Brown */
538c34d3f5SBarry Smith PetscErrorCode PetscObjectOptionsBegin_Private(PetscOptions *PetscOptionsObject,PetscObject obj)
543194b578SJed Brown {
553194b578SJed Brown   PetscErrorCode ierr;
563194b578SJed Brown   char           title[256];
573194b578SJed Brown   PetscBool      flg;
583194b578SJed Brown 
593194b578SJed Brown   PetscFunctionBegin;
603194b578SJed Brown   PetscValidHeader(obj,1);
61e55864a3SBarry Smith   PetscOptionsObject->object         = obj;
62e55864a3SBarry Smith   PetscOptionsObject->alreadyprinted = obj->optionsprinted;
63a297a907SKarl Rupp 
643194b578SJed Brown   ierr = PetscStrcmp(obj->description,obj->class_name,&flg);CHKERRQ(ierr);
653194b578SJed Brown   if (flg) {
668caf3d72SBarry Smith     ierr = PetscSNPrintf(title,sizeof(title),"%s options",obj->class_name);CHKERRQ(ierr);
673194b578SJed Brown   } else {
688caf3d72SBarry Smith     ierr = PetscSNPrintf(title,sizeof(title),"%s (%s) options",obj->description,obj->class_name);CHKERRQ(ierr);
693194b578SJed Brown   }
70e55864a3SBarry Smith   ierr = PetscOptionsBegin_Private(PetscOptionsObject,obj->comm,obj->prefix,title,obj->mansec);CHKERRQ(ierr);
713194b578SJed Brown   PetscFunctionReturn(0);
723194b578SJed Brown }
733194b578SJed Brown 
7453acd3b1SBarry Smith /*
7553acd3b1SBarry Smith      Handles adding another option to the list of options within this particular PetscOptionsBegin() PetscOptionsEnd()
7653acd3b1SBarry Smith */
7753acd3b1SBarry Smith #undef __FUNCT__
7853acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsCreate_Private"
798c34d3f5SBarry Smith static int PetscOptionsCreate_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscOptionType t,PetscOption *amsopt)
8053acd3b1SBarry Smith {
8153acd3b1SBarry Smith   int          ierr;
828c34d3f5SBarry Smith   PetscOption next;
833be6e4c3SJed Brown   PetscBool    valid;
8453acd3b1SBarry Smith 
8553acd3b1SBarry Smith   PetscFunctionBegin;
863be6e4c3SJed Brown   ierr = PetscOptionsValidKey(opt,&valid);CHKERRQ(ierr);
873be6e4c3SJed Brown   if (!valid) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"The option '%s' is not a valid key",opt);
883be6e4c3SJed Brown 
89b00a9115SJed Brown   ierr            = PetscNew(amsopt);CHKERRQ(ierr);
9053acd3b1SBarry Smith   (*amsopt)->next = 0;
9153acd3b1SBarry Smith   (*amsopt)->set  = PETSC_FALSE;
926356e834SBarry Smith   (*amsopt)->type = t;
9353acd3b1SBarry Smith   (*amsopt)->data = 0;
9461b37b28SSatish Balay 
9553acd3b1SBarry Smith   ierr = PetscStrallocpy(text,&(*amsopt)->text);CHKERRQ(ierr);
9653acd3b1SBarry Smith   ierr = PetscStrallocpy(opt,&(*amsopt)->option);CHKERRQ(ierr);
976356e834SBarry Smith   ierr = PetscStrallocpy(man,&(*amsopt)->man);CHKERRQ(ierr);
9853acd3b1SBarry Smith 
99e55864a3SBarry Smith   if (!PetscOptionsObject->next) PetscOptionsObject->next = *amsopt;
100a297a907SKarl Rupp   else {
101e55864a3SBarry Smith     next = PetscOptionsObject->next;
10253acd3b1SBarry Smith     while (next->next) next = next->next;
10353acd3b1SBarry Smith     next->next = *amsopt;
10453acd3b1SBarry Smith   }
10553acd3b1SBarry Smith   PetscFunctionReturn(0);
10653acd3b1SBarry Smith }
10753acd3b1SBarry Smith 
10853acd3b1SBarry Smith #undef __FUNCT__
109aee2cecaSBarry Smith #define __FUNCT__ "PetscScanString"
110aee2cecaSBarry Smith /*
1113fc1eb6aSBarry Smith     PetscScanString -  Gets user input via stdin from process and broadcasts to all processes
1123fc1eb6aSBarry Smith 
1133fc1eb6aSBarry Smith     Collective on MPI_Comm
1143fc1eb6aSBarry Smith 
1153fc1eb6aSBarry Smith    Input Parameters:
1163fc1eb6aSBarry Smith +     commm - communicator for the broadcast, must be PETSC_COMM_WORLD
1173fc1eb6aSBarry Smith .     n - length of the string, must be the same on all processes
1183fc1eb6aSBarry Smith -     str - location to store input
119aee2cecaSBarry Smith 
120aee2cecaSBarry Smith     Bugs:
121aee2cecaSBarry Smith .   Assumes process 0 of the given communicator has access to stdin
122aee2cecaSBarry Smith 
123aee2cecaSBarry Smith */
1243fc1eb6aSBarry Smith static PetscErrorCode PetscScanString(MPI_Comm comm,size_t n,char str[])
125aee2cecaSBarry Smith {
126330cf3c9SBarry Smith   size_t         i;
127aee2cecaSBarry Smith   char           c;
1283fc1eb6aSBarry Smith   PetscMPIInt    rank,nm;
129aee2cecaSBarry Smith   PetscErrorCode ierr;
130aee2cecaSBarry Smith 
131aee2cecaSBarry Smith   PetscFunctionBegin;
132aee2cecaSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
133aee2cecaSBarry Smith   if (!rank) {
134aee2cecaSBarry Smith     c = (char) getchar();
135aee2cecaSBarry Smith     i = 0;
136aee2cecaSBarry Smith     while (c != '\n' && i < n-1) {
137aee2cecaSBarry Smith       str[i++] = c;
138aee2cecaSBarry Smith       c = (char)getchar();
139aee2cecaSBarry Smith     }
140aee2cecaSBarry Smith     str[i] = 0;
141aee2cecaSBarry Smith   }
1424dc2109aSBarry Smith   ierr = PetscMPIIntCast(n,&nm);CHKERRQ(ierr);
1433fc1eb6aSBarry Smith   ierr = MPI_Bcast(str,nm,MPI_CHAR,0,comm);CHKERRQ(ierr);
144aee2cecaSBarry Smith   PetscFunctionReturn(0);
145aee2cecaSBarry Smith }
146aee2cecaSBarry Smith 
147ead66b60SBarry Smith #undef __FUNCT__
148ead66b60SBarry Smith #define __FUNCT__ "PetscStrdup"
1495b02f95dSBarry Smith /*
1505b02f95dSBarry Smith     This is needed because certain strings may be freed by SAWs, hence we cannot use PetscStrallocpy()
1515b02f95dSBarry Smith */
1525b02f95dSBarry Smith static PetscErrorCode  PetscStrdup(const char s[],char *t[])
1535b02f95dSBarry Smith {
1545b02f95dSBarry Smith   PetscErrorCode ierr;
1555b02f95dSBarry Smith   size_t         len;
1565b02f95dSBarry Smith   char           *tmp = 0;
1575b02f95dSBarry Smith 
1585b02f95dSBarry Smith   PetscFunctionBegin;
1595b02f95dSBarry Smith   if (s) {
1605b02f95dSBarry Smith     ierr = PetscStrlen(s,&len);CHKERRQ(ierr);
161ee48884fSBarry Smith     tmp = (char*) malloc((len+1)*sizeof(char*));
1625b02f95dSBarry Smith     if (!tmp) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"No memory to duplicate string");
1635b02f95dSBarry Smith     ierr = PetscStrcpy(tmp,s);CHKERRQ(ierr);
1645b02f95dSBarry Smith   }
1655b02f95dSBarry Smith   *t = tmp;
1665b02f95dSBarry Smith   PetscFunctionReturn(0);
1675b02f95dSBarry Smith }
1685b02f95dSBarry Smith 
1695b02f95dSBarry Smith 
170aee2cecaSBarry Smith #undef __FUNCT__
171aee2cecaSBarry Smith #define __FUNCT__ "PetscOptionsGetFromTextInput"
172aee2cecaSBarry Smith /*
1733cc1e11dSBarry Smith     PetscOptionsGetFromTextInput - Presents all the PETSc Options processed by the program so the user may change them at runtime
174aee2cecaSBarry Smith 
175aee2cecaSBarry Smith     Notes: this isn't really practical, it is just to demonstrate the principle
176aee2cecaSBarry Smith 
1777781c08eSBarry Smith     A carriage return indicates no change from the default; but this like -ksp_monitor <stdout>  the default is actually not stdout the default
1787781c08eSBarry Smith     is to do nothing so to get it to use stdout you need to type stdout. This is kind of bug?
1797781c08eSBarry Smith 
180aee2cecaSBarry Smith     Bugs:
1817781c08eSBarry Smith +    All processes must traverse through the exact same set of option queries due to the call to PetscScanString()
1823cc1e11dSBarry Smith .    Internal strings have arbitrary length and string copies are not checked that they fit into string space
183aee2cecaSBarry Smith -    Only works for PetscInt == int, PetscReal == double etc
184aee2cecaSBarry Smith 
1853cc1e11dSBarry Smith     Developer Notes: Normally the GUI that presents the options the user and retrieves the values would be running in a different
1863cc1e11dSBarry Smith      address space and communicating with the PETSc program
1873cc1e11dSBarry Smith 
188aee2cecaSBarry Smith */
1898c34d3f5SBarry Smith PetscErrorCode PetscOptionsGetFromTextInput(PetscOptions *PetscOptionsObject)
1906356e834SBarry Smith {
1916356e834SBarry Smith   PetscErrorCode ierr;
1928c34d3f5SBarry Smith   PetscOption   next = PetscOptionsObject->next;
1936356e834SBarry Smith   char           str[512];
1947781c08eSBarry Smith   PetscBool      bid;
195a4404d99SBarry Smith   PetscReal      ir,*valr;
196330cf3c9SBarry Smith   PetscInt       *vald;
197330cf3c9SBarry Smith   size_t         i;
1986356e834SBarry Smith 
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);
2971a1499c8SBarry 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);
3011a1499c8SBarry 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 
3749c1e0ce8SBarry Smith static const char *OptionsHeader = "<head>\n"
37523a1ff79SBarry Smith                                    "<script type=\"text/javascript\" src=\"http://www.mcs.anl.gov/research/projects/saws/js/jquery-1.9.1.js\"></script>\n"
37623a1ff79SBarry Smith                                    "<script type=\"text/javascript\" src=\"http://www.mcs.anl.gov/research/projects/saws/js/SAWs.js\"></script>\n"
377d1fc0251SBarry Smith                                    "<script type=\"text/javascript\" src=\"js/PETSc.js\"></script>\n"
37864bbc9efSBarry Smith                                    "<script>\n"
37964bbc9efSBarry Smith                                       "jQuery(document).ready(function() {\n"
38064bbc9efSBarry Smith                                          "PETSc.getAndDisplayDirectory(null,\"#variablesInfo\")\n"
38164bbc9efSBarry Smith                                       "})\n"
38264bbc9efSBarry Smith                                   "</script>\n"
38364bbc9efSBarry Smith                                   "</head>\n";
384426f447eSSurtai Han static const char *OptionsBodyBottom = "<div id=\"variablesInfo\" style=\"height:550px;overflow:scroll;\"></div>\n<br>\n</body>";
38564bbc9efSBarry Smith 
386d5649816SBarry Smith #undef __FUNCT__
3877781c08eSBarry Smith #define __FUNCT__ "PetscOptionsSAWsInput"
388b3506946SBarry Smith /*
3897781c08eSBarry Smith     PetscOptionsSAWsInput - Presents all the PETSc Options processed by the program so the user may change them at runtime using the SAWs
390b3506946SBarry Smith 
391b3506946SBarry Smith     Bugs:
392b3506946SBarry Smith +    All processes must traverse through the exact same set of option queries do to the call to PetscScanString()
393b3506946SBarry Smith .    Internal strings have arbitrary length and string copies are not checked that they fit into string space
394b3506946SBarry Smith -    Only works for PetscInt == int, PetscReal == double etc
395b3506946SBarry Smith 
396b3506946SBarry Smith 
397b3506946SBarry Smith */
398f7b25cf6SBarry Smith PetscErrorCode PetscOptionsSAWsInput(PetscOptions *PetscOptionsObject)
399b3506946SBarry Smith {
400b3506946SBarry Smith   PetscErrorCode ierr;
4018c34d3f5SBarry Smith   PetscOption    next     = PetscOptionsObject->next;
402d5649816SBarry Smith   static int     mancount = 0;
403b3506946SBarry Smith   char           options[16];
4047aab2a10SBarry Smith   PetscBool      changedmethod = PETSC_FALSE;
405a23eb982SSurtai Han   PetscBool      stopasking    = PETSC_FALSE;
40688a9752cSBarry Smith   char           manname[16],textname[16];
4072657e9d9SBarry Smith   char           dir[1024];
408b3506946SBarry Smith 
409b3506946SBarry Smith   /* the next line is a bug, this will only work if all processors are here, the comm passed in is ignored!!! */
410b3506946SBarry Smith   sprintf(options,"Options_%d",count++);
411a297a907SKarl Rupp 
412e55864a3SBarry Smith   PetscOptionsObject->pprefix = PetscOptionsObject->prefix; /* AMS will change this, so cannot pass prefix directly */
4131bc75a8dSBarry Smith 
414d50831c4SBarry Smith   ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s","_title");CHKERRQ(ierr);
41583355fc5SBarry Smith   PetscStackCallSAWs(SAWs_Register,(dir,&PetscOptionsObject->title,1,SAWs_READ,SAWs_STRING));
4167781c08eSBarry Smith   ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s","prefix");CHKERRQ(ierr);
41783355fc5SBarry Smith   PetscStackCallSAWs(SAWs_Register,(dir,&PetscOptionsObject->pprefix,1,SAWs_READ,SAWs_STRING));
4182657e9d9SBarry Smith   PetscStackCallSAWs(SAWs_Register,("/PETSc/Options/ChangedMethod",&changedmethod,1,SAWs_WRITE,SAWs_BOOLEAN));
419a23eb982SSurtai Han   PetscStackCallSAWs(SAWs_Register,("/PETSc/Options/StopAsking",&stopasking,1,SAWs_WRITE,SAWs_BOOLEAN));
420b3506946SBarry Smith 
421b3506946SBarry Smith   while (next) {
422d50831c4SBarry Smith     sprintf(manname,"_man_%d",mancount);
4232657e9d9SBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",manname);CHKERRQ(ierr);
4247781c08eSBarry Smith     PetscStackCallSAWs(SAWs_Register,(dir,&next->man,1,SAWs_READ,SAWs_STRING));
425d50831c4SBarry Smith     sprintf(textname,"_text_%d",mancount++);
4267781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",textname);CHKERRQ(ierr);
4277781c08eSBarry Smith     PetscStackCallSAWs(SAWs_Register,(dir,&next->text,1,SAWs_READ,SAWs_STRING));
4289f32e415SBarry Smith 
429b3506946SBarry Smith     switch (next->type) {
430b3506946SBarry Smith     case OPTION_HEAD:
431b3506946SBarry Smith       break;
432b3506946SBarry Smith     case OPTION_INT_ARRAY:
4337781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4342657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_INT));
435b3506946SBarry Smith       break;
436b3506946SBarry Smith     case OPTION_REAL_ARRAY:
4377781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4382657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_DOUBLE));
439b3506946SBarry Smith       break;
440b3506946SBarry Smith     case OPTION_INT:
4417781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4422657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_INT));
443b3506946SBarry Smith       break;
444b3506946SBarry Smith     case OPTION_REAL:
4457781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4462657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_DOUBLE));
447b3506946SBarry Smith       break;
4487781c08eSBarry Smith     case OPTION_BOOL:
4497781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4502657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_BOOLEAN));
4511ae3d29cSBarry Smith       break;
4527781c08eSBarry Smith     case OPTION_BOOL_ARRAY:
4537781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4542657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_BOOLEAN));
45571f08665SBarry Smith       break;
456b3506946SBarry Smith     case OPTION_STRING:
4577781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4587781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
4591ae3d29cSBarry Smith       break;
4601ae3d29cSBarry Smith     case OPTION_STRING_ARRAY:
4617781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4622657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_STRING));
463b3506946SBarry Smith       break;
464a264d7a6SBarry Smith     case OPTION_FLIST:
465a264d7a6SBarry Smith       {
466a264d7a6SBarry Smith       PetscInt ntext;
4677781c08eSBarry Smith       ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4687781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
469a264d7a6SBarry Smith       ierr = PetscFunctionListGet(next->flist,(const char***)&next->edata,&ntext);CHKERRQ(ierr);
470a264d7a6SBarry Smith       PetscStackCallSAWs(SAWs_Set_Legal_Variable_Values,(dir,ntext,next->edata));
471a264d7a6SBarry Smith       }
472a264d7a6SBarry Smith       break;
4731ae3d29cSBarry Smith     case OPTION_ELIST:
474a264d7a6SBarry Smith       {
475a264d7a6SBarry Smith       PetscInt ntext = next->nlist;
4767781c08eSBarry Smith       ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4777781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
478ead66b60SBarry Smith       ierr = PetscMalloc1((ntext+1),(char***)&next->edata);CHKERRQ(ierr);
4791ae3d29cSBarry Smith       ierr = PetscMemcpy(next->edata,next->list,ntext*sizeof(char*));CHKERRQ(ierr);
480a264d7a6SBarry Smith       PetscStackCallSAWs(SAWs_Set_Legal_Variable_Values,(dir,ntext,next->edata));
481a264d7a6SBarry Smith       }
482a264d7a6SBarry Smith       break;
483b3506946SBarry Smith     default:
484b3506946SBarry Smith       break;
485b3506946SBarry Smith     }
486b3506946SBarry Smith     next = next->next;
487b3506946SBarry Smith   }
488b3506946SBarry Smith 
489b3506946SBarry Smith   /* wait until accessor has unlocked the memory */
49064bbc9efSBarry Smith   PetscStackCallSAWs(SAWs_Push_Header,("index.html",OptionsHeader));
49164bbc9efSBarry Smith   PetscStackCallSAWs(SAWs_Push_Body,("index.html",2,OptionsBodyBottom));
4927aab2a10SBarry Smith   ierr = PetscSAWsBlock();CHKERRQ(ierr);
49364bbc9efSBarry Smith   PetscStackCallSAWs(SAWs_Pop_Header,("index.html"));
49464bbc9efSBarry Smith   PetscStackCallSAWs(SAWs_Pop_Body,("index.html",2));
495b3506946SBarry Smith 
49688a9752cSBarry Smith   /* determine if any values have been set in GUI */
49783355fc5SBarry Smith   next = PetscOptionsObject->next;
49888a9752cSBarry Smith   while (next) {
49988a9752cSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
500f7b25cf6SBarry Smith     PetscStackCallSAWs(SAWs_Selected,(dir,(int*)&next->set));
50188a9752cSBarry Smith     next = next->next;
50288a9752cSBarry Smith   }
50388a9752cSBarry Smith 
504b3506946SBarry Smith   /* reset counter to -2; this updates the screen with the new options for the selected method */
505f7b25cf6SBarry Smith   if (changedmethod) PetscOptionsObject->count = -2;
506b3506946SBarry Smith 
507a23eb982SSurtai Han   if (stopasking) {
508a23eb982SSurtai Han     PetscOptionsPublish      = PETSC_FALSE;
509*12655325SBarry Smith     PetscOptionsObject->count = 0;//do not ask for same thing again
510a23eb982SSurtai Han   }
511a23eb982SSurtai Han 
5129a492a5cSBarry Smith   PetscStackCallSAWs(SAWs_Delete,("/PETSc/Options"));
513b3506946SBarry Smith   PetscFunctionReturn(0);
514b3506946SBarry Smith }
515b3506946SBarry Smith #endif
516b3506946SBarry Smith 
5176356e834SBarry Smith #undef __FUNCT__
51853acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEnd_Private"
5198c34d3f5SBarry Smith PetscErrorCode PetscOptionsEnd_Private(PetscOptions *PetscOptionsObject)
52053acd3b1SBarry Smith {
52153acd3b1SBarry Smith   PetscErrorCode ierr;
5228c34d3f5SBarry Smith   PetscOption   last;
5236356e834SBarry Smith   char           option[256],value[1024],tmp[32];
524330cf3c9SBarry Smith   size_t         j;
52553acd3b1SBarry Smith 
52653acd3b1SBarry Smith   PetscFunctionBegin;
52783355fc5SBarry Smith   if (PetscOptionsObject->next) {
52883355fc5SBarry Smith     if (!PetscOptionsObject->count) {
529a264d7a6SBarry Smith #if defined(PETSC_HAVE_SAWS)
530f7b25cf6SBarry Smith       ierr = PetscOptionsSAWsInput(PetscOptionsObject);CHKERRQ(ierr);
531b3506946SBarry Smith #else
532e55864a3SBarry Smith       ierr = PetscOptionsGetFromTextInput(PetscOptionsObject);CHKERRQ(ierr);
533b3506946SBarry Smith #endif
534aee2cecaSBarry Smith     }
535aee2cecaSBarry Smith   }
5366356e834SBarry Smith 
537e55864a3SBarry Smith   ierr = PetscFree(PetscOptionsObject->title);CHKERRQ(ierr);
538e55864a3SBarry Smith   ierr = PetscFree(PetscOptionsObject->prefix);CHKERRQ(ierr);
5396356e834SBarry Smith 
540e26ddf31SBarry Smith   /* reset counter to -2; this updates the screen with the new options for the selected method */
541e55864a3SBarry Smith   if (PetscOptionsObject->changedmethod) PetscOptionsObject->count = -2;
5427a72a596SBarry Smith   /* reset alreadyprinted flag */
543e55864a3SBarry Smith   PetscOptionsObject->alreadyprinted = PETSC_FALSE;
544e55864a3SBarry Smith   if (PetscOptionsObject->object) PetscOptionsObject->object->optionsprinted = PETSC_TRUE;
545e55864a3SBarry Smith   PetscOptionsObject->object = NULL;
54653acd3b1SBarry Smith 
547e55864a3SBarry Smith   while (PetscOptionsObject->next) {
548e55864a3SBarry Smith     if (PetscOptionsObject->next->set) {
549e55864a3SBarry Smith       if (PetscOptionsObject->prefix) {
55053acd3b1SBarry Smith         ierr = PetscStrcpy(option,"-");CHKERRQ(ierr);
551e55864a3SBarry Smith         ierr = PetscStrcat(option,PetscOptionsObject->prefix);CHKERRQ(ierr);
552e55864a3SBarry Smith         ierr = PetscStrcat(option,PetscOptionsObject->next->option+1);CHKERRQ(ierr);
5536356e834SBarry Smith       } else {
554e55864a3SBarry Smith         ierr = PetscStrcpy(option,PetscOptionsObject->next->option);CHKERRQ(ierr);
5556356e834SBarry Smith       }
5566356e834SBarry Smith 
557e55864a3SBarry Smith       switch (PetscOptionsObject->next->type) {
5586356e834SBarry Smith       case OPTION_HEAD:
5596356e834SBarry Smith         break;
5606356e834SBarry Smith       case OPTION_INT_ARRAY:
561e55864a3SBarry Smith         sprintf(value,"%d",(int)((PetscInt*)PetscOptionsObject->next->data)[0]);
562e55864a3SBarry Smith         for (j=1; j<PetscOptionsObject->next->arraylength; j++) {
563e55864a3SBarry Smith           sprintf(tmp,"%d",(int)((PetscInt*)PetscOptionsObject->next->data)[j]);
5646356e834SBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5656356e834SBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5666356e834SBarry Smith         }
5676356e834SBarry Smith         break;
5686356e834SBarry Smith       case OPTION_INT:
569e55864a3SBarry Smith         sprintf(value,"%d",(int) *(PetscInt*)PetscOptionsObject->next->data);
5706356e834SBarry Smith         break;
5716356e834SBarry Smith       case OPTION_REAL:
572e55864a3SBarry Smith         sprintf(value,"%g",(double) *(PetscReal*)PetscOptionsObject->next->data);
5736356e834SBarry Smith         break;
5746356e834SBarry Smith       case OPTION_REAL_ARRAY:
575e55864a3SBarry Smith         sprintf(value,"%g",(double)((PetscReal*)PetscOptionsObject->next->data)[0]);
576e55864a3SBarry Smith         for (j=1; j<PetscOptionsObject->next->arraylength; j++) {
577e55864a3SBarry Smith           sprintf(tmp,"%g",(double)((PetscReal*)PetscOptionsObject->next->data)[j]);
5786356e834SBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5796356e834SBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5806356e834SBarry Smith         }
5816356e834SBarry Smith         break;
5827781c08eSBarry Smith       case OPTION_BOOL:
583e55864a3SBarry Smith         sprintf(value,"%d",*(int*)PetscOptionsObject->next->data);
5846356e834SBarry Smith         break;
5857781c08eSBarry Smith       case OPTION_BOOL_ARRAY:
586e55864a3SBarry Smith         sprintf(value,"%d",(int)((PetscBool*)PetscOptionsObject->next->data)[0]);
587e55864a3SBarry Smith         for (j=1; j<PetscOptionsObject->next->arraylength; j++) {
588e55864a3SBarry Smith           sprintf(tmp,"%d",(int)((PetscBool*)PetscOptionsObject->next->data)[j]);
5891ae3d29cSBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5901ae3d29cSBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5911ae3d29cSBarry Smith         }
5921ae3d29cSBarry Smith         break;
593a264d7a6SBarry Smith       case OPTION_FLIST:
5941ae3d29cSBarry Smith       case OPTION_ELIST:
595e55864a3SBarry Smith         ierr = PetscStrcpy(value,(char*)PetscOptionsObject->next->data);CHKERRQ(ierr);
5966356e834SBarry Smith         break;
5971ae3d29cSBarry Smith       case OPTION_STRING:
598e55864a3SBarry Smith         ierr = PetscStrcpy(value,(char*)PetscOptionsObject->next->data);CHKERRQ(ierr);
5991ae3d29cSBarry Smith       case OPTION_STRING_ARRAY:
600e55864a3SBarry Smith         sprintf(value,"%s",((char**)PetscOptionsObject->next->data)[0]);
601e55864a3SBarry Smith         for (j=1; j<PetscOptionsObject->next->arraylength; j++) {
602e55864a3SBarry Smith           sprintf(tmp,"%s",((char**)PetscOptionsObject->next->data)[j]);
6031ae3d29cSBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
6041ae3d29cSBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
6051ae3d29cSBarry Smith         }
6066356e834SBarry Smith         break;
6076356e834SBarry Smith       }
6086356e834SBarry Smith       ierr = PetscOptionsSetValue(option,value);CHKERRQ(ierr);
6096356e834SBarry Smith     }
610e55864a3SBarry Smith     ierr   = PetscFree(PetscOptionsObject->next->text);CHKERRQ(ierr);
611e55864a3SBarry Smith     ierr   = PetscFree(PetscOptionsObject->next->option);CHKERRQ(ierr);
612e55864a3SBarry Smith     ierr   = PetscFree(PetscOptionsObject->next->man);CHKERRQ(ierr);
613e55864a3SBarry Smith     ierr   = PetscFree(PetscOptionsObject->next->edata);CHKERRQ(ierr);
614c979a496SBarry Smith 
61583355fc5SBarry Smith     if ((PetscOptionsObject->next->type == OPTION_STRING) || (PetscOptionsObject->next->type == OPTION_FLIST) || (PetscOptionsObject->next->type == OPTION_ELIST)){
61683355fc5SBarry Smith       free(PetscOptionsObject->next->data);
617c979a496SBarry Smith     } else {
61883355fc5SBarry Smith       ierr   = PetscFree(PetscOptionsObject->next->data);CHKERRQ(ierr);
619c979a496SBarry Smith     }
6207781c08eSBarry Smith 
62183355fc5SBarry Smith     last                    = PetscOptionsObject->next;
62283355fc5SBarry Smith     PetscOptionsObject->next = PetscOptionsObject->next->next;
6236356e834SBarry Smith     ierr                    = PetscFree(last);CHKERRQ(ierr);
6246356e834SBarry Smith   }
625e55864a3SBarry Smith   PetscOptionsObject->next = 0;
62653acd3b1SBarry Smith   PetscFunctionReturn(0);
62753acd3b1SBarry Smith }
62853acd3b1SBarry Smith 
62953acd3b1SBarry Smith #undef __FUNCT__
630e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsEnum_Private"
63153acd3b1SBarry Smith /*@C
63253acd3b1SBarry Smith    PetscOptionsEnum - Gets the enum value for a particular option in the database.
63353acd3b1SBarry Smith 
6343f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
63553acd3b1SBarry Smith 
63653acd3b1SBarry Smith    Input Parameters:
63753acd3b1SBarry Smith +  opt - option name
63853acd3b1SBarry Smith .  text - short string that describes the option
63953acd3b1SBarry Smith .  man - manual page with additional information on option
64053acd3b1SBarry Smith .  list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null
6410fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
6420fdccdaeSBarry Smith $                 PetscOptionsEnum(..., obj->value,&object->value,...) or
6430fdccdaeSBarry Smith $                 value = defaultvalue
6440fdccdaeSBarry Smith $                 PetscOptionsEnum(..., value,&value,&flg);
6450fdccdaeSBarry Smith $                 if (flg) {
64653acd3b1SBarry Smith 
64753acd3b1SBarry Smith    Output Parameter:
64853acd3b1SBarry Smith +  value - the  value to return
649b32e0204SMatthew G Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
65053acd3b1SBarry Smith 
65153acd3b1SBarry Smith    Level: beginner
65253acd3b1SBarry Smith 
65353acd3b1SBarry Smith    Concepts: options database
65453acd3b1SBarry Smith 
65553acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
65653acd3b1SBarry Smith 
65753acd3b1SBarry Smith           list is usually something like PCASMTypes or some other predefined list of enum names
65853acd3b1SBarry Smith 
65953acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
660acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
661acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
66253acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
66353acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
664acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
665a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
66653acd3b1SBarry Smith @*/
6678c34d3f5SBarry Smith PetscErrorCode  PetscOptionsEnum_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],const char *const *list,PetscEnum currentvalue,PetscEnum *value,PetscBool  *set)
66853acd3b1SBarry Smith {
66953acd3b1SBarry Smith   PetscErrorCode ierr;
67053acd3b1SBarry Smith   PetscInt       ntext = 0;
671aa5bb8c0SSatish Balay   PetscInt       tval;
672ace3abfcSBarry Smith   PetscBool      tflg;
67353acd3b1SBarry Smith 
67453acd3b1SBarry Smith   PetscFunctionBegin;
67553acd3b1SBarry Smith   while (list[ntext++]) {
676e32f2f54SBarry Smith     if (ntext > 50) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument appears to be wrong or have more than 50 entries");
67753acd3b1SBarry Smith   }
678e32f2f54SBarry Smith   if (ntext < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument must have at least two entries: typename and type prefix");
67953acd3b1SBarry Smith   ntext -= 3;
680e55864a3SBarry Smith   ierr   = PetscOptionsEList_Private(PetscOptionsObject,opt,text,man,list,ntext,list[currentvalue],&tval,&tflg);CHKERRQ(ierr);
681aa5bb8c0SSatish Balay   /* with PETSC_USE_64BIT_INDICES sizeof(PetscInt) != sizeof(PetscEnum) */
682aa5bb8c0SSatish Balay   if (tflg) *value = (PetscEnum)tval;
683aa5bb8c0SSatish Balay   if (set)  *set   = tflg;
68453acd3b1SBarry Smith   PetscFunctionReturn(0);
68553acd3b1SBarry Smith }
68653acd3b1SBarry Smith 
68753acd3b1SBarry Smith /* -------------------------------------------------------------------------------------------------------------*/
68853acd3b1SBarry Smith #undef __FUNCT__
689e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsInt_Private"
69053acd3b1SBarry Smith /*@C
69153acd3b1SBarry Smith    PetscOptionsInt - Gets the integer value for a particular option in the database.
69253acd3b1SBarry Smith 
6933f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
69453acd3b1SBarry Smith 
69553acd3b1SBarry Smith    Input Parameters:
69653acd3b1SBarry Smith +  opt - option name
69753acd3b1SBarry Smith .  text - short string that describes the option
69853acd3b1SBarry Smith .  man - manual page with additional information on option
6990fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
7000fdccdaeSBarry Smith $                 PetscOptionsInt(..., obj->value,&object->value,...) or
7010fdccdaeSBarry Smith $                 value = defaultvalue
7020fdccdaeSBarry Smith $                 PetscOptionsInt(..., value,&value,&flg);
7030fdccdaeSBarry Smith $                 if (flg) {
70453acd3b1SBarry Smith 
70553acd3b1SBarry Smith    Output Parameter:
70653acd3b1SBarry Smith +  value - the integer value to return
70753acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
70853acd3b1SBarry Smith 
70953acd3b1SBarry Smith    Level: beginner
71053acd3b1SBarry Smith 
71153acd3b1SBarry Smith    Concepts: options database^has int
71253acd3b1SBarry Smith 
71353acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
71453acd3b1SBarry Smith 
71553acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
716acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
717acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
71853acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
71953acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
720acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
721a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
72253acd3b1SBarry Smith @*/
7238c34d3f5SBarry Smith PetscErrorCode  PetscOptionsInt_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscInt currentvalue,PetscInt *value,PetscBool  *set)
72453acd3b1SBarry Smith {
72553acd3b1SBarry Smith   PetscErrorCode ierr;
7268c34d3f5SBarry Smith   PetscOption    amsopt;
727*12655325SBarry Smith   PetscBool      wasset;
72853acd3b1SBarry Smith 
72953acd3b1SBarry Smith   PetscFunctionBegin;
730e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
731e55864a3SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_INT,&amsopt);CHKERRQ(ierr);
7326356e834SBarry Smith     ierr = PetscMalloc(sizeof(PetscInt),&amsopt->data);CHKERRQ(ierr);
733*12655325SBarry Smith     *(PetscInt*)amsopt->data = currentvalue;
7343e211508SBarry Smith 
735*12655325SBarry Smith     ierr = PetscOptionsGetInt(PetscOptionsObject->prefix,opt,&currentvalue,&wasset);CHKERRQ(ierr);
7363e211508SBarry Smith     if (wasset) {
737*12655325SBarry Smith       *(PetscInt*)amsopt->data = currentvalue;
7383e211508SBarry Smith     }
739af6d86caSBarry Smith   }
740e55864a3SBarry Smith   ierr = PetscOptionsGetInt(PetscOptionsObject->prefix,opt,value,set);CHKERRQ(ierr);
741e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
7421a1499c8SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%d>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,currentvalue,text,ManSection(man));CHKERRQ(ierr);
74353acd3b1SBarry Smith   }
74453acd3b1SBarry Smith   PetscFunctionReturn(0);
74553acd3b1SBarry Smith }
74653acd3b1SBarry Smith 
74753acd3b1SBarry Smith #undef __FUNCT__
748e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsString_Private"
74953acd3b1SBarry Smith /*@C
75053acd3b1SBarry Smith    PetscOptionsString - Gets the string value for a particular option in the database.
75153acd3b1SBarry Smith 
7523f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
75353acd3b1SBarry Smith 
75453acd3b1SBarry Smith    Input Parameters:
75553acd3b1SBarry Smith +  opt - option name
75653acd3b1SBarry Smith .  text - short string that describes the option
75753acd3b1SBarry Smith .  man - manual page with additional information on option
7580fdccdaeSBarry Smith .  currentvalue - the current value; caller is responsible for setting this value correctly. This is not used to set value
759bcbf2dc5SJed Brown -  len - length of the result string including null terminator
76053acd3b1SBarry Smith 
76153acd3b1SBarry Smith    Output Parameter:
76253acd3b1SBarry Smith +  value - the value to return
76353acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
76453acd3b1SBarry Smith 
76553acd3b1SBarry Smith    Level: beginner
76653acd3b1SBarry Smith 
76753acd3b1SBarry Smith    Concepts: options database^has int
76853acd3b1SBarry Smith 
76953acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
77053acd3b1SBarry Smith 
7717fccdfe4SBarry 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).
7727fccdfe4SBarry Smith 
77353acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
774acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
775acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsReal(), PetscOptionsBool(),
77653acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
77753acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
778acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
779a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
78053acd3b1SBarry Smith @*/
7818c34d3f5SBarry Smith PetscErrorCode  PetscOptionsString_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],const char currentvalue[],char value[],size_t len,PetscBool  *set)
78253acd3b1SBarry Smith {
78353acd3b1SBarry Smith   PetscErrorCode ierr;
7848c34d3f5SBarry Smith   PetscOption   amsopt;
78553acd3b1SBarry Smith 
78653acd3b1SBarry Smith   PetscFunctionBegin;
7871a1499c8SBarry Smith   if (!PetscOptionsObject->count) {
7881a1499c8SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
78964facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
7900fdccdaeSBarry Smith     ierr = PetscStrdup(currentvalue ? currentvalue : "",(char**)&amsopt->data);CHKERRQ(ierr);
791af6d86caSBarry Smith   }
792e55864a3SBarry Smith   ierr = PetscOptionsGetString(PetscOptionsObject->prefix,opt,value,len,set);CHKERRQ(ierr);
793e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
7941a1499c8SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%s>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,currentvalue,text,ManSection(man));CHKERRQ(ierr);
79553acd3b1SBarry Smith   }
79653acd3b1SBarry Smith   PetscFunctionReturn(0);
79753acd3b1SBarry Smith }
79853acd3b1SBarry Smith 
79953acd3b1SBarry Smith #undef __FUNCT__
800e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsReal_Private"
80153acd3b1SBarry Smith /*@C
80253acd3b1SBarry Smith    PetscOptionsReal - Gets the PetscReal value for a particular option in the database.
80353acd3b1SBarry Smith 
8043f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
80553acd3b1SBarry Smith 
80653acd3b1SBarry Smith    Input Parameters:
80753acd3b1SBarry Smith +  opt - option name
80853acd3b1SBarry Smith .  text - short string that describes the option
80953acd3b1SBarry Smith .  man - manual page with additional information on option
8100fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
8110fdccdaeSBarry Smith $                 PetscOptionsReal(..., obj->value,&object->value,...) or
8120fdccdaeSBarry Smith $                 value = defaultvalue
8130fdccdaeSBarry Smith $                 PetscOptionsReal(..., value,&value,&flg);
8140fdccdaeSBarry Smith $                 if (flg) {
81553acd3b1SBarry Smith 
81653acd3b1SBarry Smith    Output Parameter:
81753acd3b1SBarry Smith +  value - the value to return
81853acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
81953acd3b1SBarry Smith 
82053acd3b1SBarry Smith    Level: beginner
82153acd3b1SBarry Smith 
82253acd3b1SBarry Smith    Concepts: options database^has int
82353acd3b1SBarry Smith 
82453acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
82553acd3b1SBarry Smith 
82653acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
827acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
828acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
82953acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
83053acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
831acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
832a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
83353acd3b1SBarry Smith @*/
8348c34d3f5SBarry Smith PetscErrorCode  PetscOptionsReal_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscReal currentvalue,PetscReal *value,PetscBool  *set)
83553acd3b1SBarry Smith {
83653acd3b1SBarry Smith   PetscErrorCode ierr;
8378c34d3f5SBarry Smith   PetscOption   amsopt;
83853acd3b1SBarry Smith 
83953acd3b1SBarry Smith   PetscFunctionBegin;
840e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
841e55864a3SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_REAL,&amsopt);CHKERRQ(ierr);
842538aa990SBarry Smith     ierr = PetscMalloc(sizeof(PetscReal),&amsopt->data);CHKERRQ(ierr);
843a297a907SKarl Rupp 
8440fdccdaeSBarry Smith     *(PetscReal*)amsopt->data = currentvalue;
845538aa990SBarry Smith   }
8461a1499c8SBarry Smith   ierr = PetscOptionsGetReal(PetscOptionsObject->prefix,opt,value,set);CHKERRQ(ierr);
8471a1499c8SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
8481a1499c8SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%g>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,(double)currentvalue,text,ManSection(man));CHKERRQ(ierr);
84953acd3b1SBarry Smith   }
85053acd3b1SBarry Smith   PetscFunctionReturn(0);
85153acd3b1SBarry Smith }
85253acd3b1SBarry Smith 
85353acd3b1SBarry Smith #undef __FUNCT__
854e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsScalar_Private"
85553acd3b1SBarry Smith /*@C
85653acd3b1SBarry Smith    PetscOptionsScalar - Gets the scalar value for a particular option in the database.
85753acd3b1SBarry Smith 
8583f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
85953acd3b1SBarry Smith 
86053acd3b1SBarry Smith    Input Parameters:
86153acd3b1SBarry Smith +  opt - option name
86253acd3b1SBarry Smith .  text - short string that describes the option
86353acd3b1SBarry Smith .  man - manual page with additional information on option
8640fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
8650fdccdaeSBarry Smith $                 PetscOptionsScalar(..., obj->value,&object->value,...) or
8660fdccdaeSBarry Smith $                 value = defaultvalue
8670fdccdaeSBarry Smith $                 PetscOptionsScalar(..., value,&value,&flg);
8680fdccdaeSBarry Smith $                 if (flg) {
8690fdccdaeSBarry Smith 
87053acd3b1SBarry Smith 
87153acd3b1SBarry Smith    Output Parameter:
87253acd3b1SBarry Smith +  value - the value to return
87353acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
87453acd3b1SBarry Smith 
87553acd3b1SBarry Smith    Level: beginner
87653acd3b1SBarry Smith 
87753acd3b1SBarry Smith    Concepts: options database^has int
87853acd3b1SBarry Smith 
87953acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
88053acd3b1SBarry Smith 
88153acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
882acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
883acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
88453acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
88553acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
886acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
887a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
88853acd3b1SBarry Smith @*/
8898c34d3f5SBarry Smith PetscErrorCode  PetscOptionsScalar_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscScalar currentvalue,PetscScalar *value,PetscBool  *set)
89053acd3b1SBarry Smith {
89153acd3b1SBarry Smith   PetscErrorCode ierr;
89253acd3b1SBarry Smith 
89353acd3b1SBarry Smith   PetscFunctionBegin;
89453acd3b1SBarry Smith #if !defined(PETSC_USE_COMPLEX)
8950fdccdaeSBarry Smith   ierr = PetscOptionsReal(opt,text,man,currentvalue,value,set);CHKERRQ(ierr);
89653acd3b1SBarry Smith #else
897e55864a3SBarry Smith   ierr = PetscOptionsGetScalar(PetscOptionsObject->prefix,opt,value,set);CHKERRQ(ierr);
89853acd3b1SBarry Smith #endif
89953acd3b1SBarry Smith   PetscFunctionReturn(0);
90053acd3b1SBarry Smith }
90153acd3b1SBarry Smith 
90253acd3b1SBarry Smith #undef __FUNCT__
903e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsName_Private"
90453acd3b1SBarry Smith /*@C
90590d69ab7SBarry 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
90690d69ab7SBarry Smith                       its value is set to false.
90753acd3b1SBarry Smith 
9083f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
90953acd3b1SBarry Smith 
91053acd3b1SBarry Smith    Input Parameters:
91153acd3b1SBarry Smith +  opt - option name
91253acd3b1SBarry Smith .  text - short string that describes the option
91353acd3b1SBarry Smith -  man - manual page with additional information on option
91453acd3b1SBarry Smith 
91553acd3b1SBarry Smith    Output Parameter:
91653acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
91753acd3b1SBarry Smith 
91853acd3b1SBarry Smith    Level: beginner
91953acd3b1SBarry Smith 
92053acd3b1SBarry Smith    Concepts: options database^has int
92153acd3b1SBarry Smith 
92253acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
92353acd3b1SBarry Smith 
92453acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
925acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
926acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
92753acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
92853acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
929acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
930a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
93153acd3b1SBarry Smith @*/
9328c34d3f5SBarry Smith PetscErrorCode  PetscOptionsName_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool  *flg)
93353acd3b1SBarry Smith {
93453acd3b1SBarry Smith   PetscErrorCode ierr;
9358c34d3f5SBarry Smith   PetscOption   amsopt;
93653acd3b1SBarry Smith 
93753acd3b1SBarry Smith   PetscFunctionBegin;
938e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
939e55864a3SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
940ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
941a297a907SKarl Rupp 
942ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
9431ae3d29cSBarry Smith   }
944e55864a3SBarry Smith   ierr = PetscOptionsHasName(PetscOptionsObject->prefix,opt,flg);CHKERRQ(ierr);
945e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
946e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
94753acd3b1SBarry Smith   }
94853acd3b1SBarry Smith   PetscFunctionReturn(0);
94953acd3b1SBarry Smith }
95053acd3b1SBarry Smith 
95153acd3b1SBarry Smith #undef __FUNCT__
952e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsFList_Private"
95353acd3b1SBarry Smith /*@C
954a264d7a6SBarry Smith      PetscOptionsFList - Puts a list of option values that a single one may be selected from
95553acd3b1SBarry Smith 
9563f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
95753acd3b1SBarry Smith 
95853acd3b1SBarry Smith    Input Parameters:
95953acd3b1SBarry Smith +  opt - option name
96053acd3b1SBarry Smith .  text - short string that describes the option
96153acd3b1SBarry Smith .  man - manual page with additional information on option
96253acd3b1SBarry Smith .  list - the possible choices
9630fdccdaeSBarry Smith .  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with
9640fdccdaeSBarry Smith $                 PetscOptionsFlist(..., obj->value,value,len,&flg);
9650fdccdaeSBarry Smith $                 if (flg) {
9663cc1e11dSBarry Smith -  len - the length of the character array value
96753acd3b1SBarry Smith 
96853acd3b1SBarry Smith    Output Parameter:
96953acd3b1SBarry Smith +  value - the value to return
97053acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
97153acd3b1SBarry Smith 
97253acd3b1SBarry Smith    Level: intermediate
97353acd3b1SBarry Smith 
97453acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
97553acd3b1SBarry Smith 
97653acd3b1SBarry Smith    See PetscOptionsEList() for when the choices are given in a string array
97753acd3b1SBarry Smith 
97853acd3b1SBarry Smith    To get a listing of all currently specified options,
97988c29154SBarry Smith     see PetscOptionsView() or PetscOptionsGetAll()
98053acd3b1SBarry Smith 
981eabe10d7SBarry Smith    Developer Note: This cannot check for invalid selection because of things like MATAIJ that are not included in the list
982eabe10d7SBarry Smith 
98353acd3b1SBarry Smith    Concepts: options database^list
98453acd3b1SBarry Smith 
98553acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
986acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
98753acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
98853acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
989acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
990a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsEnum()
99153acd3b1SBarry Smith @*/
9928c34d3f5SBarry Smith PetscErrorCode  PetscOptionsFList_Private(PetscOptions *PetscOptionsObject,const char opt[],const char ltext[],const char man[],PetscFunctionList list,const char currentvalue[],char value[],size_t len,PetscBool  *set)
99353acd3b1SBarry Smith {
99453acd3b1SBarry Smith   PetscErrorCode ierr;
9958c34d3f5SBarry Smith   PetscOption   amsopt;
99653acd3b1SBarry Smith 
99753acd3b1SBarry Smith   PetscFunctionBegin;
9981a1499c8SBarry Smith   if (!PetscOptionsObject->count) {
9991a1499c8SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,ltext,man,OPTION_FLIST,&amsopt);CHKERRQ(ierr);
100064facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
10010fdccdaeSBarry Smith     ierr = PetscStrdup(currentvalue ? currentvalue : "",(char**)&amsopt->data);CHKERRQ(ierr);
10023cc1e11dSBarry Smith     amsopt->flist = list;
10033cc1e11dSBarry Smith   }
10041a1499c8SBarry Smith   ierr = PetscOptionsGetString(PetscOptionsObject->prefix,opt,value,len,set);CHKERRQ(ierr);
10051a1499c8SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
10061a1499c8SBarry Smith     ierr = PetscFunctionListPrintTypes(PetscOptionsObject->comm,stdout,PetscOptionsObject->prefix,opt,ltext,man,list,currentvalue);CHKERRQ(ierr);CHKERRQ(ierr);
100753acd3b1SBarry Smith   }
100853acd3b1SBarry Smith   PetscFunctionReturn(0);
100953acd3b1SBarry Smith }
101053acd3b1SBarry Smith 
101153acd3b1SBarry Smith #undef __FUNCT__
1012e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsEList_Private"
101353acd3b1SBarry Smith /*@C
101453acd3b1SBarry Smith      PetscOptionsEList - Puts a list of option values that a single one may be selected from
101553acd3b1SBarry Smith 
10163f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
101753acd3b1SBarry Smith 
101853acd3b1SBarry Smith    Input Parameters:
101953acd3b1SBarry Smith +  opt - option name
102053acd3b1SBarry Smith .  ltext - short string that describes the option
102153acd3b1SBarry Smith .  man - manual page with additional information on option
1022a264d7a6SBarry Smith .  list - the possible choices (one of these must be selected, anything else is invalid)
102353acd3b1SBarry Smith .  ntext - number of choices
10240fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with
10250fdccdaeSBarry Smith $                 PetscOptionsElist(..., obj->value,&value,&flg);
10260fdccdaeSBarry Smith $                 if (flg) {
10270fdccdaeSBarry Smith 
102853acd3b1SBarry Smith 
102953acd3b1SBarry Smith    Output Parameter:
103053acd3b1SBarry Smith +  value - the index of the value to return
103153acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
103253acd3b1SBarry Smith 
103353acd3b1SBarry Smith    Level: intermediate
103453acd3b1SBarry Smith 
103553acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
103653acd3b1SBarry Smith 
1037a264d7a6SBarry Smith    See PetscOptionsFList() for when the choices are given in a PetscFunctionList()
103853acd3b1SBarry Smith 
103953acd3b1SBarry Smith    Concepts: options database^list
104053acd3b1SBarry Smith 
104153acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1042acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
104353acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
104453acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1045acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1046a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEnum()
104753acd3b1SBarry Smith @*/
10488c34d3f5SBarry Smith PetscErrorCode  PetscOptionsEList_Private(PetscOptions *PetscOptionsObject,const char opt[],const char ltext[],const char man[],const char *const *list,PetscInt ntext,const char currentvalue[],PetscInt *value,PetscBool  *set)
104953acd3b1SBarry Smith {
105053acd3b1SBarry Smith   PetscErrorCode ierr;
105153acd3b1SBarry Smith   PetscInt       i;
10528c34d3f5SBarry Smith   PetscOption   amsopt;
105353acd3b1SBarry Smith 
105453acd3b1SBarry Smith   PetscFunctionBegin;
10551a1499c8SBarry Smith   if (!PetscOptionsObject->count) {
10561a1499c8SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,ltext,man,OPTION_ELIST,&amsopt);CHKERRQ(ierr);
105764facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
10580fdccdaeSBarry Smith     ierr = PetscStrdup(currentvalue ? currentvalue : "",(char**)&amsopt->data);CHKERRQ(ierr);
10591ae3d29cSBarry Smith     amsopt->list  = list;
10601ae3d29cSBarry Smith     amsopt->nlist = ntext;
10611ae3d29cSBarry Smith   }
10621a1499c8SBarry Smith   ierr = PetscOptionsGetEList(PetscOptionsObject->prefix,opt,list,ntext,value,set);CHKERRQ(ierr);
10631a1499c8SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
10641a1499c8SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%s> (choose one of)",PetscOptionsObject->prefix?PetscOptionsObject->prefix:"",opt+1,currentvalue);CHKERRQ(ierr);
106553acd3b1SBarry Smith     for (i=0; i<ntext; i++) {
1066e55864a3SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," %s",list[i]);CHKERRQ(ierr);
106753acd3b1SBarry Smith     }
1068e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," (%s)\n",ManSection(man));CHKERRQ(ierr);
106953acd3b1SBarry Smith   }
107053acd3b1SBarry Smith   PetscFunctionReturn(0);
107153acd3b1SBarry Smith }
107253acd3b1SBarry Smith 
107353acd3b1SBarry Smith #undef __FUNCT__
1074e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsBoolGroupBegin_Private"
107553acd3b1SBarry Smith /*@C
1076acfcf0e5SJed Brown      PetscOptionsBoolGroupBegin - First in a series of logical queries on the options database for
1077d5649816SBarry Smith        which at most a single value can be true.
107853acd3b1SBarry Smith 
10793f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
108053acd3b1SBarry Smith 
108153acd3b1SBarry Smith    Input Parameters:
108253acd3b1SBarry Smith +  opt - option name
108353acd3b1SBarry Smith .  text - short string that describes the option
108453acd3b1SBarry Smith -  man - manual page with additional information on option
108553acd3b1SBarry Smith 
108653acd3b1SBarry Smith    Output Parameter:
108753acd3b1SBarry Smith .  flg - whether that option was set or not
108853acd3b1SBarry Smith 
108953acd3b1SBarry Smith    Level: intermediate
109053acd3b1SBarry Smith 
109153acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
109253acd3b1SBarry Smith 
1093acfcf0e5SJed Brown    Must be followed by 0 or more PetscOptionsBoolGroup()s and PetscOptionsBoolGroupEnd()
109453acd3b1SBarry Smith 
109553acd3b1SBarry Smith     Concepts: options database^logical group
109653acd3b1SBarry Smith 
109753acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1098acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
109953acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
110053acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1101acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1102a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
110353acd3b1SBarry Smith @*/
11048c34d3f5SBarry Smith PetscErrorCode  PetscOptionsBoolGroupBegin_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool  *flg)
110553acd3b1SBarry Smith {
110653acd3b1SBarry Smith   PetscErrorCode ierr;
11078c34d3f5SBarry Smith   PetscOption   amsopt;
110853acd3b1SBarry Smith 
110953acd3b1SBarry Smith   PetscFunctionBegin;
1110e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
111183355fc5SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1112ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1113a297a907SKarl Rupp 
1114ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
11151ae3d29cSBarry Smith   }
111668b16fdaSBarry Smith   *flg = PETSC_FALSE;
1117e55864a3SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject->prefix,opt,flg,NULL);CHKERRQ(ierr);
1118e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1119e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  Pick at most one of -------------\n");CHKERRQ(ierr);
1120e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"    -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
112153acd3b1SBarry Smith   }
112253acd3b1SBarry Smith   PetscFunctionReturn(0);
112353acd3b1SBarry Smith }
112453acd3b1SBarry Smith 
112553acd3b1SBarry Smith #undef __FUNCT__
1126e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsBoolGroup_Private"
112753acd3b1SBarry Smith /*@C
1128acfcf0e5SJed Brown      PetscOptionsBoolGroup - One in a series of logical queries on the options database for
1129d5649816SBarry Smith        which at most a single value can be true.
113053acd3b1SBarry Smith 
11313f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
113253acd3b1SBarry Smith 
113353acd3b1SBarry Smith    Input Parameters:
113453acd3b1SBarry Smith +  opt - option name
113553acd3b1SBarry Smith .  text - short string that describes the option
113653acd3b1SBarry Smith -  man - manual page with additional information on option
113753acd3b1SBarry Smith 
113853acd3b1SBarry Smith    Output Parameter:
113953acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
114053acd3b1SBarry Smith 
114153acd3b1SBarry Smith    Level: intermediate
114253acd3b1SBarry Smith 
114353acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
114453acd3b1SBarry Smith 
1145acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin() and preceded a PetscOptionsBoolGroupEnd()
114653acd3b1SBarry Smith 
114753acd3b1SBarry Smith     Concepts: options database^logical group
114853acd3b1SBarry Smith 
114953acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1150acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
115153acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
115253acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1153acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1154a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
115553acd3b1SBarry Smith @*/
11568c34d3f5SBarry Smith PetscErrorCode  PetscOptionsBoolGroup_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool  *flg)
115753acd3b1SBarry Smith {
115853acd3b1SBarry Smith   PetscErrorCode ierr;
11598c34d3f5SBarry Smith   PetscOption   amsopt;
116053acd3b1SBarry Smith 
116153acd3b1SBarry Smith   PetscFunctionBegin;
1162e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
116383355fc5SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1164ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1165a297a907SKarl Rupp 
1166ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
11671ae3d29cSBarry Smith   }
116817326d04SJed Brown   *flg = PETSC_FALSE;
1169e55864a3SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject->prefix,opt,flg,NULL);CHKERRQ(ierr);
1170e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1171e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"    -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
117253acd3b1SBarry Smith   }
117353acd3b1SBarry Smith   PetscFunctionReturn(0);
117453acd3b1SBarry Smith }
117553acd3b1SBarry Smith 
117653acd3b1SBarry Smith #undef __FUNCT__
1177e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsBoolGroupEnd_Private"
117853acd3b1SBarry Smith /*@C
1179acfcf0e5SJed Brown      PetscOptionsBoolGroupEnd - Last in a series of logical queries on the options database for
1180d5649816SBarry Smith        which at most a single value can be true.
118153acd3b1SBarry Smith 
11823f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
118353acd3b1SBarry Smith 
118453acd3b1SBarry Smith    Input Parameters:
118553acd3b1SBarry Smith +  opt - option name
118653acd3b1SBarry Smith .  text - short string that describes the option
118753acd3b1SBarry Smith -  man - manual page with additional information on option
118853acd3b1SBarry Smith 
118953acd3b1SBarry Smith    Output Parameter:
119053acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
119153acd3b1SBarry Smith 
119253acd3b1SBarry Smith    Level: intermediate
119353acd3b1SBarry Smith 
119453acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
119553acd3b1SBarry Smith 
1196acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin()
119753acd3b1SBarry Smith 
119853acd3b1SBarry Smith     Concepts: options database^logical group
119953acd3b1SBarry Smith 
120053acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1201acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
120253acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
120353acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1204acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1205a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
120653acd3b1SBarry Smith @*/
12078c34d3f5SBarry Smith PetscErrorCode  PetscOptionsBoolGroupEnd_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool  *flg)
120853acd3b1SBarry Smith {
120953acd3b1SBarry Smith   PetscErrorCode ierr;
12108c34d3f5SBarry Smith   PetscOption   amsopt;
121153acd3b1SBarry Smith 
121253acd3b1SBarry Smith   PetscFunctionBegin;
1213e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
121483355fc5SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1215ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1216a297a907SKarl Rupp 
1217ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
12181ae3d29cSBarry Smith   }
121917326d04SJed Brown   *flg = PETSC_FALSE;
1220e55864a3SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject->prefix,opt,flg,NULL);CHKERRQ(ierr);
1221e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1222e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"    -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
122353acd3b1SBarry Smith   }
122453acd3b1SBarry Smith   PetscFunctionReturn(0);
122553acd3b1SBarry Smith }
122653acd3b1SBarry Smith 
122753acd3b1SBarry Smith #undef __FUNCT__
1228e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsBool_Private"
122953acd3b1SBarry Smith /*@C
1230acfcf0e5SJed Brown    PetscOptionsBool - Determines if a particular option is in the database with a true or false
123153acd3b1SBarry Smith 
12323f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
123353acd3b1SBarry Smith 
123453acd3b1SBarry Smith    Input Parameters:
123553acd3b1SBarry Smith +  opt - option name
123653acd3b1SBarry Smith .  text - short string that describes the option
1237868c398cSBarry Smith .  man - manual page with additional information on option
123894ae4db5SBarry Smith -  currentvalue - the current value
123953acd3b1SBarry Smith 
124053acd3b1SBarry Smith    Output Parameter:
124153acd3b1SBarry Smith .  flg - PETSC_TRUE or PETSC_FALSE
124253acd3b1SBarry Smith .  set - PETSC_TRUE if found, else PETSC_FALSE
124353acd3b1SBarry Smith 
124453acd3b1SBarry Smith    Level: beginner
124553acd3b1SBarry Smith 
124653acd3b1SBarry Smith    Concepts: options database^logical
124753acd3b1SBarry Smith 
124853acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
124953acd3b1SBarry Smith 
125053acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
1251acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
1252acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
125353acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
125453acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1255acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1256a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
125753acd3b1SBarry Smith @*/
12588c34d3f5SBarry Smith PetscErrorCode  PetscOptionsBool_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool currentvalue,PetscBool  *flg,PetscBool  *set)
125953acd3b1SBarry Smith {
126053acd3b1SBarry Smith   PetscErrorCode ierr;
1261ace3abfcSBarry Smith   PetscBool      iset;
12628c34d3f5SBarry Smith   PetscOption   amsopt;
126353acd3b1SBarry Smith 
126453acd3b1SBarry Smith   PetscFunctionBegin;
1265e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
1266e55864a3SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1267ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1268a297a907SKarl Rupp 
126994ae4db5SBarry Smith     *(PetscBool*)amsopt->data = currentvalue;
1270af6d86caSBarry Smith   }
12711a1499c8SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject->prefix,opt,flg,&iset);CHKERRQ(ierr);
127253acd3b1SBarry Smith   if (set) *set = iset;
12731a1499c8SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
127494ae4db5SBarry Smith     const char *v = PetscBools[currentvalue];
12751a1499c8SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s: <%s> %s (%s)\n",PetscOptionsObject->prefix?PetscOptionsObject->prefix:"",opt+1,v,text,ManSection(man));CHKERRQ(ierr);
127653acd3b1SBarry Smith   }
127753acd3b1SBarry Smith   PetscFunctionReturn(0);
127853acd3b1SBarry Smith }
127953acd3b1SBarry Smith 
128053acd3b1SBarry Smith #undef __FUNCT__
1281e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsRealArray_Private"
128253acd3b1SBarry Smith /*@C
128353acd3b1SBarry Smith    PetscOptionsRealArray - Gets an array of double values for a particular
128453acd3b1SBarry Smith    option in the database. The values must be separated with commas with
128553acd3b1SBarry Smith    no intervening spaces.
128653acd3b1SBarry Smith 
12873f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
128853acd3b1SBarry Smith 
128953acd3b1SBarry Smith    Input Parameters:
129053acd3b1SBarry Smith +  opt - the option one is seeking
129153acd3b1SBarry Smith .  text - short string describing option
129253acd3b1SBarry Smith .  man - manual page for option
129353acd3b1SBarry Smith -  nmax - maximum number of values
129453acd3b1SBarry Smith 
129553acd3b1SBarry Smith    Output Parameter:
129653acd3b1SBarry Smith +  value - location to copy values
129753acd3b1SBarry Smith .  nmax - actual number of values found
129853acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
129953acd3b1SBarry Smith 
130053acd3b1SBarry Smith    Level: beginner
130153acd3b1SBarry Smith 
130253acd3b1SBarry Smith    Notes:
130353acd3b1SBarry Smith    The user should pass in an array of doubles
130453acd3b1SBarry Smith 
130553acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
130653acd3b1SBarry Smith 
130753acd3b1SBarry Smith    Concepts: options database^array of strings
130853acd3b1SBarry Smith 
130953acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1310acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
131153acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
131253acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1313acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1314a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
131553acd3b1SBarry Smith @*/
13168c34d3f5SBarry Smith PetscErrorCode  PetscOptionsRealArray_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscReal value[],PetscInt *n,PetscBool  *set)
131753acd3b1SBarry Smith {
131853acd3b1SBarry Smith   PetscErrorCode ierr;
131953acd3b1SBarry Smith   PetscInt       i;
13208c34d3f5SBarry Smith   PetscOption   amsopt;
132153acd3b1SBarry Smith 
132253acd3b1SBarry Smith   PetscFunctionBegin;
1323e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
1324e26ddf31SBarry Smith     PetscReal *vals;
1325e26ddf31SBarry Smith 
1326e55864a3SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_REAL_ARRAY,&amsopt);CHKERRQ(ierr);
1327e55864a3SBarry Smith     ierr = PetscMalloc((*n)*sizeof(PetscReal),&amsopt->data);CHKERRQ(ierr);
1328e26ddf31SBarry Smith     vals = (PetscReal*)amsopt->data;
1329e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1330e26ddf31SBarry Smith     amsopt->arraylength = *n;
1331e26ddf31SBarry Smith   }
1332e55864a3SBarry Smith   ierr = PetscOptionsGetRealArray(PetscOptionsObject->prefix,opt,value,n,set);CHKERRQ(ierr);
1333e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1334a519f713SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%g",PetscOptionsObject->prefix?PetscOptionsObject->prefix:"",opt+1,(double)value[0]);CHKERRQ(ierr);
133553acd3b1SBarry Smith     for (i=1; i<*n; i++) {
1336e55864a3SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%g",(double)value[i]);CHKERRQ(ierr);
133753acd3b1SBarry Smith     }
1338e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
133953acd3b1SBarry Smith   }
134053acd3b1SBarry Smith   PetscFunctionReturn(0);
134153acd3b1SBarry Smith }
134253acd3b1SBarry Smith 
134353acd3b1SBarry Smith 
134453acd3b1SBarry Smith #undef __FUNCT__
1345e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsIntArray_Private"
134653acd3b1SBarry Smith /*@C
134753acd3b1SBarry Smith    PetscOptionsIntArray - Gets an array of integers for a particular
1348b32a342fSShri Abhyankar    option in the database.
134953acd3b1SBarry Smith 
13503f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
135153acd3b1SBarry Smith 
135253acd3b1SBarry Smith    Input Parameters:
135353acd3b1SBarry Smith +  opt - the option one is seeking
135453acd3b1SBarry Smith .  text - short string describing option
135553acd3b1SBarry Smith .  man - manual page for option
1356f8a50e2bSBarry Smith -  n - maximum number of values
135753acd3b1SBarry Smith 
135853acd3b1SBarry Smith    Output Parameter:
135953acd3b1SBarry Smith +  value - location to copy values
1360f8a50e2bSBarry Smith .  n - actual number of values found
136153acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
136253acd3b1SBarry Smith 
136353acd3b1SBarry Smith    Level: beginner
136453acd3b1SBarry Smith 
136553acd3b1SBarry Smith    Notes:
1366b32a342fSShri Abhyankar    The array can be passed as
1367b32a342fSShri Abhyankar    a comma seperated list:                                 0,1,2,3,4,5,6,7
13680fd488f5SShri Abhyankar    a range (start-end+1):                                  0-8
13690fd488f5SShri Abhyankar    a range with given increment (start-end+1:inc):         0-7:2
13700fd488f5SShri Abhyankar    a combination of values and ranges seperated by commas: 0,1-8,8-15:2
1371b32a342fSShri Abhyankar 
1372b32a342fSShri Abhyankar    There must be no intervening spaces between the values.
137353acd3b1SBarry Smith 
137453acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
137553acd3b1SBarry Smith 
1376b32a342fSShri Abhyankar    Concepts: options database^array of ints
137753acd3b1SBarry Smith 
137853acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1379acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
138053acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
138153acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1382acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1383a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsRealArray()
138453acd3b1SBarry Smith @*/
13858c34d3f5SBarry Smith PetscErrorCode  PetscOptionsIntArray_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscInt value[],PetscInt *n,PetscBool  *set)
138653acd3b1SBarry Smith {
138753acd3b1SBarry Smith   PetscErrorCode ierr;
138853acd3b1SBarry Smith   PetscInt       i;
13898c34d3f5SBarry Smith   PetscOption   amsopt;
139053acd3b1SBarry Smith 
139153acd3b1SBarry Smith   PetscFunctionBegin;
1392e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
1393e26ddf31SBarry Smith     PetscInt *vals;
1394e26ddf31SBarry Smith 
1395e55864a3SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_INT_ARRAY,&amsopt);CHKERRQ(ierr);
1396854ce69bSBarry Smith     ierr = PetscMalloc1(*n,(PetscInt**)&amsopt->data);CHKERRQ(ierr);
1397e26ddf31SBarry Smith     vals = (PetscInt*)amsopt->data;
1398e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1399e26ddf31SBarry Smith     amsopt->arraylength = *n;
1400e26ddf31SBarry Smith   }
1401e55864a3SBarry Smith   ierr = PetscOptionsGetIntArray(PetscOptionsObject->prefix,opt,value,n,set);CHKERRQ(ierr);
1402e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1403e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%d",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,value[0]);CHKERRQ(ierr);
140453acd3b1SBarry Smith     for (i=1; i<*n; i++) {
1405e55864a3SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%d",value[i]);CHKERRQ(ierr);
140653acd3b1SBarry Smith     }
1407e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
140853acd3b1SBarry Smith   }
140953acd3b1SBarry Smith   PetscFunctionReturn(0);
141053acd3b1SBarry Smith }
141153acd3b1SBarry Smith 
141253acd3b1SBarry Smith #undef __FUNCT__
1413e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsStringArray_Private"
141453acd3b1SBarry Smith /*@C
141553acd3b1SBarry Smith    PetscOptionsStringArray - Gets an array of string values for a particular
141653acd3b1SBarry Smith    option in the database. The values must be separated with commas with
141753acd3b1SBarry Smith    no intervening spaces.
141853acd3b1SBarry Smith 
14193f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
142053acd3b1SBarry Smith 
142153acd3b1SBarry Smith    Input Parameters:
142253acd3b1SBarry Smith +  opt - the option one is seeking
142353acd3b1SBarry Smith .  text - short string describing option
142453acd3b1SBarry Smith .  man - manual page for option
142553acd3b1SBarry Smith -  nmax - maximum number of strings
142653acd3b1SBarry Smith 
142753acd3b1SBarry Smith    Output Parameter:
142853acd3b1SBarry Smith +  value - location to copy strings
142953acd3b1SBarry Smith .  nmax - actual number of strings found
143053acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
143153acd3b1SBarry Smith 
143253acd3b1SBarry Smith    Level: beginner
143353acd3b1SBarry Smith 
143453acd3b1SBarry Smith    Notes:
143553acd3b1SBarry Smith    The user should pass in an array of pointers to char, to hold all the
143653acd3b1SBarry Smith    strings returned by this function.
143753acd3b1SBarry Smith 
143853acd3b1SBarry Smith    The user is responsible for deallocating the strings that are
143953acd3b1SBarry Smith    returned. The Fortran interface for this routine is not supported.
144053acd3b1SBarry Smith 
144153acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
144253acd3b1SBarry Smith 
144353acd3b1SBarry Smith    Concepts: options database^array of strings
144453acd3b1SBarry Smith 
144553acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1446acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
144753acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
144853acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1449acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1450a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
145153acd3b1SBarry Smith @*/
14528c34d3f5SBarry Smith PetscErrorCode  PetscOptionsStringArray_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],char *value[],PetscInt *nmax,PetscBool  *set)
145353acd3b1SBarry Smith {
145453acd3b1SBarry Smith   PetscErrorCode ierr;
14558c34d3f5SBarry Smith   PetscOption   amsopt;
145653acd3b1SBarry Smith 
145753acd3b1SBarry Smith   PetscFunctionBegin;
1458e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
1459e55864a3SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_STRING_ARRAY,&amsopt);CHKERRQ(ierr);
1460854ce69bSBarry Smith     ierr = PetscMalloc1(*nmax,(char**)&amsopt->data);CHKERRQ(ierr);
1461a297a907SKarl Rupp 
14621ae3d29cSBarry Smith     amsopt->arraylength = *nmax;
14631ae3d29cSBarry Smith   }
1464e55864a3SBarry Smith   ierr = PetscOptionsGetStringArray(PetscOptionsObject->prefix,opt,value,nmax,set);CHKERRQ(ierr);
1465e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1466e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <string1,string2,...>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
146753acd3b1SBarry Smith   }
146853acd3b1SBarry Smith   PetscFunctionReturn(0);
146953acd3b1SBarry Smith }
147053acd3b1SBarry Smith 
1471e2446a98SMatthew Knepley #undef __FUNCT__
1472e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsBoolArray_Private"
1473e2446a98SMatthew Knepley /*@C
1474acfcf0e5SJed Brown    PetscOptionsBoolArray - Gets an array of logical values (true or false) for a particular
1475e2446a98SMatthew Knepley    option in the database. The values must be separated with commas with
1476e2446a98SMatthew Knepley    no intervening spaces.
1477e2446a98SMatthew Knepley 
14783f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
1479e2446a98SMatthew Knepley 
1480e2446a98SMatthew Knepley    Input Parameters:
1481e2446a98SMatthew Knepley +  opt - the option one is seeking
1482e2446a98SMatthew Knepley .  text - short string describing option
1483e2446a98SMatthew Knepley .  man - manual page for option
1484e2446a98SMatthew Knepley -  nmax - maximum number of values
1485e2446a98SMatthew Knepley 
1486e2446a98SMatthew Knepley    Output Parameter:
1487e2446a98SMatthew Knepley +  value - location to copy values
1488e2446a98SMatthew Knepley .  nmax - actual number of values found
1489e2446a98SMatthew Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
1490e2446a98SMatthew Knepley 
1491e2446a98SMatthew Knepley    Level: beginner
1492e2446a98SMatthew Knepley 
1493e2446a98SMatthew Knepley    Notes:
1494e2446a98SMatthew Knepley    The user should pass in an array of doubles
1495e2446a98SMatthew Knepley 
1496e2446a98SMatthew Knepley    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
1497e2446a98SMatthew Knepley 
1498e2446a98SMatthew Knepley    Concepts: options database^array of strings
1499e2446a98SMatthew Knepley 
1500e2446a98SMatthew Knepley .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1501acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
1502e2446a98SMatthew Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1503e2446a98SMatthew Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1504acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1505a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
1506e2446a98SMatthew Knepley @*/
15078c34d3f5SBarry Smith PetscErrorCode  PetscOptionsBoolArray_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool value[],PetscInt *n,PetscBool *set)
1508e2446a98SMatthew Knepley {
1509e2446a98SMatthew Knepley   PetscErrorCode ierr;
1510e2446a98SMatthew Knepley   PetscInt       i;
15118c34d3f5SBarry Smith   PetscOption   amsopt;
1512e2446a98SMatthew Knepley 
1513e2446a98SMatthew Knepley   PetscFunctionBegin;
1514e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
1515ace3abfcSBarry Smith     PetscBool *vals;
15161ae3d29cSBarry Smith 
151783355fc5SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL_ARRAY,&amsopt);CHKERRQ(ierr);
15181a1499c8SBarry Smith     ierr = PetscMalloc1(*n,(PetscBool**)&amsopt->data);CHKERRQ(ierr);
1519ace3abfcSBarry Smith     vals = (PetscBool*)amsopt->data;
15201ae3d29cSBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
15211ae3d29cSBarry Smith     amsopt->arraylength = *n;
15221ae3d29cSBarry Smith   }
1523e55864a3SBarry Smith   ierr = PetscOptionsGetBoolArray(PetscOptionsObject->prefix,opt,value,n,set);CHKERRQ(ierr);
1524e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1525e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%d",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,value[0]);CHKERRQ(ierr);
1526e2446a98SMatthew Knepley     for (i=1; i<*n; i++) {
1527e55864a3SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%d",value[i]);CHKERRQ(ierr);
1528e2446a98SMatthew Knepley     }
1529e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
1530e2446a98SMatthew Knepley   }
1531e2446a98SMatthew Knepley   PetscFunctionReturn(0);
1532e2446a98SMatthew Knepley }
1533e2446a98SMatthew Knepley 
15348cc676e6SMatthew G Knepley #undef __FUNCT__
1535e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsViewer_Private"
15368cc676e6SMatthew G Knepley /*@C
1537d1da0b69SBarry Smith    PetscOptionsViewer - Gets a viewer appropriate for the type indicated by the user
15388cc676e6SMatthew G Knepley 
15398cc676e6SMatthew G Knepley    Logically Collective on the communicator passed in PetscOptionsBegin()
15408cc676e6SMatthew G Knepley 
15418cc676e6SMatthew G Knepley    Input Parameters:
15428cc676e6SMatthew G Knepley +  opt - option name
15438cc676e6SMatthew G Knepley .  text - short string that describes the option
15448cc676e6SMatthew G Knepley -  man - manual page with additional information on option
15458cc676e6SMatthew G Knepley 
15468cc676e6SMatthew G Knepley    Output Parameter:
15478cc676e6SMatthew G Knepley +  viewer - the viewer
15488cc676e6SMatthew G Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
15498cc676e6SMatthew G Knepley 
15508cc676e6SMatthew G Knepley    Level: beginner
15518cc676e6SMatthew G Knepley 
15528cc676e6SMatthew G Knepley    Concepts: options database^has int
15538cc676e6SMatthew G Knepley 
15548cc676e6SMatthew G Knepley    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
15558cc676e6SMatthew G Knepley 
15565a7113b9SPatrick Sanan    See PetscOptionsGetViewer() for the format of the supplied viewer and its options
15578cc676e6SMatthew G Knepley 
15588cc676e6SMatthew G Knepley .seealso: PetscOptionsGetViewer(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
15598cc676e6SMatthew G Knepley           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
15608cc676e6SMatthew G Knepley           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
15618cc676e6SMatthew G Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
15628cc676e6SMatthew G Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
15638cc676e6SMatthew G Knepley           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1564a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
15658cc676e6SMatthew G Knepley @*/
15668c34d3f5SBarry Smith PetscErrorCode  PetscOptionsViewer_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscViewer *viewer,PetscViewerFormat *format,PetscBool  *set)
15678cc676e6SMatthew G Knepley {
15688cc676e6SMatthew G Knepley   PetscErrorCode ierr;
15698c34d3f5SBarry Smith   PetscOption    amsopt;
15708cc676e6SMatthew G Knepley 
15718cc676e6SMatthew G Knepley   PetscFunctionBegin;
15721a1499c8SBarry Smith   if (!PetscOptionsObject->count) {
15731a1499c8SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
157464facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
15755b02f95dSBarry Smith     ierr = PetscStrdup("",(char**)&amsopt->data);CHKERRQ(ierr);
15768cc676e6SMatthew G Knepley   }
1577e55864a3SBarry Smith   ierr = PetscOptionsGetViewer(PetscOptionsObject->comm,PetscOptionsObject->prefix,opt,viewer,format,set);CHKERRQ(ierr);
1578e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1579e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%s>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,"",text,ManSection(man));CHKERRQ(ierr);
15808cc676e6SMatthew G Knepley   }
15818cc676e6SMatthew G Knepley   PetscFunctionReturn(0);
15828cc676e6SMatthew G Knepley }
15838cc676e6SMatthew G Knepley 
158453acd3b1SBarry Smith 
158553acd3b1SBarry Smith #undef __FUNCT__
158653acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsHead"
158753acd3b1SBarry Smith /*@C
1588b52f573bSBarry Smith      PetscOptionsHead - Puts a heading before listing any more published options. Used, for example,
158953acd3b1SBarry Smith             in KSPSetFromOptions_GMRES().
159053acd3b1SBarry Smith 
15913f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
159253acd3b1SBarry Smith 
159353acd3b1SBarry Smith    Input Parameter:
159453acd3b1SBarry Smith .   head - the heading text
159553acd3b1SBarry Smith 
159653acd3b1SBarry Smith 
159753acd3b1SBarry Smith    Level: intermediate
159853acd3b1SBarry Smith 
159953acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
160053acd3b1SBarry Smith 
1601b52f573bSBarry Smith           Can be followed by a call to PetscOptionsTail() in the same function.
160253acd3b1SBarry Smith 
160353acd3b1SBarry Smith    Concepts: options database^subheading
160453acd3b1SBarry Smith 
160553acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1606acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
160753acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
160853acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1609acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1610a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
161153acd3b1SBarry Smith @*/
16128c34d3f5SBarry Smith PetscErrorCode  PetscOptionsHead(PetscOptions *PetscOptionsObject,const char head[])
161353acd3b1SBarry Smith {
161453acd3b1SBarry Smith   PetscErrorCode ierr;
161553acd3b1SBarry Smith 
161653acd3b1SBarry Smith   PetscFunctionBegin;
1617e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1618e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  %s\n",head);CHKERRQ(ierr);
161953acd3b1SBarry Smith   }
162053acd3b1SBarry Smith   PetscFunctionReturn(0);
162153acd3b1SBarry Smith }
162253acd3b1SBarry Smith 
162353acd3b1SBarry Smith 
162453acd3b1SBarry Smith 
162553acd3b1SBarry Smith 
162653acd3b1SBarry Smith 
162753acd3b1SBarry Smith 
1628