xref: /petsc/src/sys/objects/aoptions.c (revision 5a856986583887c326abe5dfd149e8184a29cd80)
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 
9af0996ceSBarry Smith #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 /*
2353acd3b1SBarry Smith     Handles setting up the data structure in a call to PetscOptionsBegin()
2453acd3b1SBarry Smith */
254416b707SBarry Smith PetscErrorCode PetscOptionsBegin_Private(PetscOptionItems *PetscOptionsObject,MPI_Comm comm,const char prefix[],const char title[],const char mansec[])
2653acd3b1SBarry Smith {
2753acd3b1SBarry Smith   PetscErrorCode ierr;
2853acd3b1SBarry Smith 
2953acd3b1SBarry Smith   PetscFunctionBegin;
300eb63584SBarry Smith   if (!PetscOptionsObject->alreadyprinted) {
319de0f6ecSBarry Smith     if (!PetscOptionsHelpPrintedSingleton) {
329de0f6ecSBarry Smith       ierr = PetscOptionsHelpPrintedCreate(&PetscOptionsHelpPrintedSingleton);CHKERRQ(ierr);
339de0f6ecSBarry Smith     }
349de0f6ecSBarry Smith     ierr = PetscOptionsHelpPrintedCheck(PetscOptionsHelpPrintedSingleton,prefix,title,&PetscOptionsObject->alreadyprinted);CHKERRQ(ierr);
350eb63584SBarry Smith   }
36e55864a3SBarry Smith   PetscOptionsObject->next          = 0;
37e55864a3SBarry Smith   PetscOptionsObject->comm          = comm;
38e55864a3SBarry Smith   PetscOptionsObject->changedmethod = PETSC_FALSE;
39a297a907SKarl Rupp 
40e55864a3SBarry Smith   ierr = PetscStrallocpy(prefix,&PetscOptionsObject->prefix);CHKERRQ(ierr);
41e55864a3SBarry Smith   ierr = PetscStrallocpy(title,&PetscOptionsObject->title);CHKERRQ(ierr);
4253acd3b1SBarry Smith 
432d747510SLisandro Dalcin   ierr = PetscOptionsHasHelp(PetscOptionsObject->options,&PetscOptionsObject->printhelp);CHKERRQ(ierr);
44e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1) {
45e55864a3SBarry Smith     if (!PetscOptionsObject->alreadyprinted) {
4653acd3b1SBarry Smith       ierr = (*PetscHelpPrintf)(comm,"%s -------------------------------------------------\n",title);CHKERRQ(ierr);
4753acd3b1SBarry Smith     }
4861b37b28SSatish Balay   }
4953acd3b1SBarry Smith   PetscFunctionReturn(0);
5053acd3b1SBarry Smith }
5153acd3b1SBarry Smith 
523194b578SJed Brown /*
533194b578SJed Brown     Handles setting up the data structure in a call to PetscObjectOptionsBegin()
543194b578SJed Brown */
554416b707SBarry Smith PetscErrorCode PetscObjectOptionsBegin_Private(PetscOptionItems *PetscOptionsObject,PetscObject obj)
563194b578SJed Brown {
573194b578SJed Brown   PetscErrorCode ierr;
583194b578SJed Brown   char           title[256];
593194b578SJed Brown   PetscBool      flg;
603194b578SJed Brown 
613194b578SJed Brown   PetscFunctionBegin;
623194b578SJed Brown   PetscValidHeader(obj,1);
63e55864a3SBarry Smith   PetscOptionsObject->object         = obj;
64e55864a3SBarry Smith   PetscOptionsObject->alreadyprinted = obj->optionsprinted;
65a297a907SKarl Rupp 
663194b578SJed Brown   ierr = PetscStrcmp(obj->description,obj->class_name,&flg);CHKERRQ(ierr);
673194b578SJed Brown   if (flg) {
688caf3d72SBarry Smith     ierr = PetscSNPrintf(title,sizeof(title),"%s options",obj->class_name);CHKERRQ(ierr);
693194b578SJed Brown   } else {
708caf3d72SBarry Smith     ierr = PetscSNPrintf(title,sizeof(title),"%s (%s) options",obj->description,obj->class_name);CHKERRQ(ierr);
713194b578SJed Brown   }
72e55864a3SBarry Smith   ierr = PetscOptionsBegin_Private(PetscOptionsObject,obj->comm,obj->prefix,title,obj->mansec);CHKERRQ(ierr);
733194b578SJed Brown   PetscFunctionReturn(0);
743194b578SJed Brown }
753194b578SJed Brown 
7653acd3b1SBarry Smith /*
7753acd3b1SBarry Smith      Handles adding another option to the list of options within this particular PetscOptionsBegin() PetscOptionsEnd()
7853acd3b1SBarry Smith */
794416b707SBarry Smith static int PetscOptionItemCreate_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscOptionType t,PetscOptionItem *amsopt)
8053acd3b1SBarry Smith {
8153acd3b1SBarry Smith   int             ierr;
824416b707SBarry Smith   PetscOptionItem 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 
108aee2cecaSBarry Smith /*
1093fc1eb6aSBarry Smith     PetscScanString -  Gets user input via stdin from process and broadcasts to all processes
1103fc1eb6aSBarry Smith 
111d083f849SBarry Smith     Collective
1123fc1eb6aSBarry Smith 
1133fc1eb6aSBarry Smith    Input Parameters:
1143fc1eb6aSBarry Smith +     commm - communicator for the broadcast, must be PETSC_COMM_WORLD
1153fc1eb6aSBarry Smith .     n - length of the string, must be the same on all processes
1163fc1eb6aSBarry Smith -     str - location to store input
117aee2cecaSBarry Smith 
118aee2cecaSBarry Smith     Bugs:
119aee2cecaSBarry Smith .   Assumes process 0 of the given communicator has access to stdin
120aee2cecaSBarry Smith 
121aee2cecaSBarry Smith */
1223fc1eb6aSBarry Smith static PetscErrorCode PetscScanString(MPI_Comm comm,size_t n,char str[])
123aee2cecaSBarry Smith {
124330cf3c9SBarry Smith   size_t         i;
125aee2cecaSBarry Smith   char           c;
1263fc1eb6aSBarry Smith   PetscMPIInt    rank,nm;
127aee2cecaSBarry Smith   PetscErrorCode ierr;
128aee2cecaSBarry Smith 
129aee2cecaSBarry Smith   PetscFunctionBegin;
130aee2cecaSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
131aee2cecaSBarry Smith   if (!rank) {
132aee2cecaSBarry Smith     c = (char) getchar();
133aee2cecaSBarry Smith     i = 0;
134aee2cecaSBarry Smith     while (c != '\n' && i < n-1) {
135aee2cecaSBarry Smith       str[i++] = c;
136aee2cecaSBarry Smith       c = (char)getchar();
137aee2cecaSBarry Smith     }
138aee2cecaSBarry Smith     str[i] = 0;
139aee2cecaSBarry Smith   }
1404dc2109aSBarry Smith   ierr = PetscMPIIntCast(n,&nm);CHKERRQ(ierr);
1413fc1eb6aSBarry Smith   ierr = MPI_Bcast(str,nm,MPI_CHAR,0,comm);CHKERRQ(ierr);
142aee2cecaSBarry Smith   PetscFunctionReturn(0);
143aee2cecaSBarry Smith }
144aee2cecaSBarry Smith 
1455b02f95dSBarry Smith /*
1465b02f95dSBarry Smith     This is needed because certain strings may be freed by SAWs, hence we cannot use PetscStrallocpy()
1475b02f95dSBarry Smith */
1485b02f95dSBarry Smith static PetscErrorCode  PetscStrdup(const char s[],char *t[])
1495b02f95dSBarry Smith {
1505b02f95dSBarry Smith   PetscErrorCode ierr;
1515b02f95dSBarry Smith   size_t         len;
1525b02f95dSBarry Smith   char           *tmp = 0;
1535b02f95dSBarry Smith 
1545b02f95dSBarry Smith   PetscFunctionBegin;
1555b02f95dSBarry Smith   if (s) {
1565b02f95dSBarry Smith     ierr = PetscStrlen(s,&len);CHKERRQ(ierr);
157f416af30SBarry Smith     tmp = (char*) malloc((len+1)*sizeof(char));
1585b02f95dSBarry Smith     if (!tmp) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"No memory to duplicate string");
1595b02f95dSBarry Smith     ierr = PetscStrcpy(tmp,s);CHKERRQ(ierr);
1605b02f95dSBarry Smith   }
1615b02f95dSBarry Smith   *t = tmp;
1625b02f95dSBarry Smith   PetscFunctionReturn(0);
1635b02f95dSBarry Smith }
1645b02f95dSBarry Smith 
1655b02f95dSBarry Smith 
166aee2cecaSBarry Smith /*
1673cc1e11dSBarry Smith     PetscOptionsGetFromTextInput - Presents all the PETSc Options processed by the program so the user may change them at runtime
168aee2cecaSBarry Smith 
16995452b02SPatrick Sanan     Notes:
17095452b02SPatrick Sanan     this isn't really practical, it is just to demonstrate the principle
171aee2cecaSBarry Smith 
1727781c08eSBarry Smith     A carriage return indicates no change from the default; but this like -ksp_monitor <stdout>  the default is actually not stdout the default
1737781c08eSBarry Smith     is to do nothing so to get it to use stdout you need to type stdout. This is kind of bug?
1747781c08eSBarry Smith 
175aee2cecaSBarry Smith     Bugs:
1767781c08eSBarry Smith +    All processes must traverse through the exact same set of option queries due to the call to PetscScanString()
1773cc1e11dSBarry Smith .    Internal strings have arbitrary length and string copies are not checked that they fit into string space
178aee2cecaSBarry Smith -    Only works for PetscInt == int, PetscReal == double etc
179aee2cecaSBarry Smith 
18095452b02SPatrick Sanan     Developer Notes:
18195452b02SPatrick Sanan     Normally the GUI that presents the options the user and retrieves the values would be running in a different
1823cc1e11dSBarry Smith      address space and communicating with the PETSc program
1833cc1e11dSBarry Smith 
184aee2cecaSBarry Smith */
1854416b707SBarry Smith PetscErrorCode PetscOptionsGetFromTextInput(PetscOptionItems *PetscOptionsObject)
1866356e834SBarry Smith {
1876356e834SBarry Smith   PetscErrorCode  ierr;
1884416b707SBarry Smith   PetscOptionItem next = PetscOptionsObject->next;
1896356e834SBarry Smith   char            str[512];
1907781c08eSBarry Smith   PetscBool       bid;
191a4404d99SBarry Smith   PetscReal       ir,*valr;
192330cf3c9SBarry Smith   PetscInt        *vald;
193330cf3c9SBarry Smith   size_t          i;
1946356e834SBarry Smith 
1952a409bb0SBarry Smith   PetscFunctionBegin;
196e55864a3SBarry Smith   ierr = (*PetscPrintf)(PETSC_COMM_WORLD,"%s -------------------------------------------------\n",PetscOptionsObject->title);CHKERRQ(ierr);
1976356e834SBarry Smith   while (next) {
1986356e834SBarry Smith     switch (next->type) {
1996356e834SBarry Smith     case OPTION_HEAD:
2006356e834SBarry Smith       break;
201e26ddf31SBarry Smith     case OPTION_INT_ARRAY:
202e55864a3SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",next->option+1);CHKERRQ(ierr);
203e26ddf31SBarry Smith       vald = (PetscInt*) next->data;
204e26ddf31SBarry Smith       for (i=0; i<next->arraylength; i++) {
205e26ddf31SBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"%d",vald[i]);CHKERRQ(ierr);
206e26ddf31SBarry Smith         if (i < next->arraylength-1) {
207e26ddf31SBarry Smith           ierr = PetscPrintf(PETSC_COMM_WORLD,",");CHKERRQ(ierr);
208e26ddf31SBarry Smith         }
209e26ddf31SBarry Smith       }
210e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,">: %s (%s) ",next->text,next->man);CHKERRQ(ierr);
211e26ddf31SBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
212e26ddf31SBarry Smith       if (str[0]) {
213e26ddf31SBarry Smith         PetscToken token;
214e26ddf31SBarry Smith         PetscInt   n=0,nmax = next->arraylength,*dvalue = (PetscInt*)next->data,start,end;
215e26ddf31SBarry Smith         size_t     len;
216e26ddf31SBarry Smith         char       *value;
217ace3abfcSBarry Smith         PetscBool  foundrange;
218e26ddf31SBarry Smith 
219e26ddf31SBarry Smith         next->set = PETSC_TRUE;
220e26ddf31SBarry Smith         value     = str;
221e26ddf31SBarry Smith         ierr      = PetscTokenCreate(value,',',&token);CHKERRQ(ierr);
222e26ddf31SBarry Smith         ierr      = PetscTokenFind(token,&value);CHKERRQ(ierr);
223e26ddf31SBarry Smith         while (n < nmax) {
224e26ddf31SBarry Smith           if (!value) break;
225e26ddf31SBarry Smith 
226e26ddf31SBarry Smith           /* look for form  d-D where d and D are integers */
227e26ddf31SBarry Smith           foundrange = PETSC_FALSE;
228e26ddf31SBarry Smith           ierr       = PetscStrlen(value,&len);CHKERRQ(ierr);
229e26ddf31SBarry Smith           if (value[0] == '-') i=2;
230e26ddf31SBarry Smith           else i=1;
231330cf3c9SBarry Smith           for (;i<len; i++) {
232e26ddf31SBarry Smith             if (value[i] == '-') {
233e32f2f54SBarry Smith               if (i == len-1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Error in %D-th array entry %s\n",n,value);
234e26ddf31SBarry Smith               value[i] = 0;
235cfbddea1SSatish Balay               ierr     = PetscOptionsStringToInt(value,&start);CHKERRQ(ierr);
236cfbddea1SSatish Balay               ierr     = PetscOptionsStringToInt(value+i+1,&end);CHKERRQ(ierr);
237e32f2f54SBarry 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);
238e32f2f54SBarry 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);
239e26ddf31SBarry Smith               for (; start<end; start++) {
240e26ddf31SBarry Smith                 *dvalue = start; dvalue++;n++;
241e26ddf31SBarry Smith               }
242e26ddf31SBarry Smith               foundrange = PETSC_TRUE;
243e26ddf31SBarry Smith               break;
244e26ddf31SBarry Smith             }
245e26ddf31SBarry Smith           }
246e26ddf31SBarry Smith           if (!foundrange) {
247cfbddea1SSatish Balay             ierr = PetscOptionsStringToInt(value,dvalue);CHKERRQ(ierr);
248e26ddf31SBarry Smith             dvalue++;
249e26ddf31SBarry Smith             n++;
250e26ddf31SBarry Smith           }
251e26ddf31SBarry Smith           ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
252e26ddf31SBarry Smith         }
2538c74ee41SBarry Smith         ierr = PetscTokenDestroy(&token);CHKERRQ(ierr);
254e26ddf31SBarry Smith       }
255e26ddf31SBarry Smith       break;
256e26ddf31SBarry Smith     case OPTION_REAL_ARRAY:
257e55864a3SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",next->option+1);CHKERRQ(ierr);
258e26ddf31SBarry Smith       valr = (PetscReal*) next->data;
259e26ddf31SBarry Smith       for (i=0; i<next->arraylength; i++) {
260e26ddf31SBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"%g",valr[i]);CHKERRQ(ierr);
261e26ddf31SBarry Smith         if (i < next->arraylength-1) {
262e26ddf31SBarry Smith           ierr = PetscPrintf(PETSC_COMM_WORLD,",");CHKERRQ(ierr);
263e26ddf31SBarry Smith         }
264e26ddf31SBarry Smith       }
265e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,">: %s (%s) ",next->text,next->man);CHKERRQ(ierr);
266e26ddf31SBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
267e26ddf31SBarry Smith       if (str[0]) {
268e26ddf31SBarry Smith         PetscToken token;
269e26ddf31SBarry Smith         PetscInt   n = 0,nmax = next->arraylength;
270e26ddf31SBarry Smith         PetscReal  *dvalue = (PetscReal*)next->data;
271e26ddf31SBarry Smith         char       *value;
272e26ddf31SBarry Smith 
273e26ddf31SBarry Smith         next->set = PETSC_TRUE;
274e26ddf31SBarry Smith         value     = str;
275e26ddf31SBarry Smith         ierr      = PetscTokenCreate(value,',',&token);CHKERRQ(ierr);
276e26ddf31SBarry Smith         ierr      = PetscTokenFind(token,&value);CHKERRQ(ierr);
277e26ddf31SBarry Smith         while (n < nmax) {
278e26ddf31SBarry Smith           if (!value) break;
279cfbddea1SSatish Balay           ierr = PetscOptionsStringToReal(value,dvalue);CHKERRQ(ierr);
280e26ddf31SBarry Smith           dvalue++;
281e26ddf31SBarry Smith           n++;
282e26ddf31SBarry Smith           ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
283e26ddf31SBarry Smith         }
2848c74ee41SBarry Smith         ierr = PetscTokenDestroy(&token);CHKERRQ(ierr);
285e26ddf31SBarry Smith       }
286e26ddf31SBarry Smith       break;
2876356e834SBarry Smith     case OPTION_INT:
288e55864a3SBarry 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);
2893fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
2903fc1eb6aSBarry Smith       if (str[0]) {
291d25d7f95SJed Brown #if defined(PETSC_SIZEOF_LONG_LONG)
292d25d7f95SJed Brown         long long lid;
293d25d7f95SJed Brown         sscanf(str,"%lld",&lid);
2941a1499c8SBarry 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);
295c272547aSJed Brown #else
296d25d7f95SJed Brown         long  lid;
297d25d7f95SJed Brown         sscanf(str,"%ld",&lid);
2981a1499c8SBarry 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);
299c272547aSJed Brown #endif
300a297a907SKarl Rupp 
301d25d7f95SJed Brown         next->set = PETSC_TRUE;
302d25d7f95SJed Brown         *((PetscInt*)next->data) = (PetscInt)lid;
303aee2cecaSBarry Smith       }
304aee2cecaSBarry Smith       break;
305aee2cecaSBarry Smith     case OPTION_REAL:
306e55864a3SBarry 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);
3073fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3083fc1eb6aSBarry Smith       if (str[0]) {
309ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
310a4404d99SBarry Smith         sscanf(str,"%e",&ir);
311570b7f6dSBarry Smith #elif defined(PETSC_USE_REAL___FP16)
312570b7f6dSBarry Smith         float irtemp;
313570b7f6dSBarry Smith         sscanf(str,"%e",&irtemp);
314570b7f6dSBarry Smith         ir = irtemp;
315ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL_DOUBLE)
316aee2cecaSBarry Smith         sscanf(str,"%le",&ir);
317ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL___FLOAT128)
318d9822059SBarry Smith         ir = strtoflt128(str,0);
319d9822059SBarry Smith #else
320513dbe71SLisandro Dalcin         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Unknown scalar type");
321a4404d99SBarry Smith #endif
322aee2cecaSBarry Smith         next->set                 = PETSC_TRUE;
323aee2cecaSBarry Smith         *((PetscReal*)next->data) = ir;
324aee2cecaSBarry Smith       }
325aee2cecaSBarry Smith       break;
3267781c08eSBarry Smith     case OPTION_BOOL:
32783355fc5SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <%s>: %s (%s) ",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",next->option+1,*(PetscBool*)next->data ? "true": "false",next->text,next->man);CHKERRQ(ierr);
3287781c08eSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3297781c08eSBarry Smith       if (str[0]) {
3307781c08eSBarry Smith         ierr = PetscOptionsStringToBool(str,&bid);CHKERRQ(ierr);
3317781c08eSBarry Smith         next->set = PETSC_TRUE;
3327781c08eSBarry Smith         *((PetscBool*)next->data) = bid;
3337781c08eSBarry Smith       }
3347781c08eSBarry Smith       break;
335aee2cecaSBarry Smith     case OPTION_STRING:
336e55864a3SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <%s>: %s (%s) ",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",next->option+1,(char*)next->data,next->text,next->man);CHKERRQ(ierr);
3373fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3383fc1eb6aSBarry Smith       if (str[0]) {
339aee2cecaSBarry Smith         next->set = PETSC_TRUE;
34064facd6cSBarry Smith         /* must use system malloc since SAWs may free this */
3415b02f95dSBarry Smith         ierr = PetscStrdup(str,(char**)&next->data);CHKERRQ(ierr);
3426356e834SBarry Smith       }
3436356e834SBarry Smith       break;
344a264d7a6SBarry Smith     case OPTION_FLIST:
34544ef3d73SBarry Smith       ierr = PetscFunctionListPrintTypes(PETSC_COMM_WORLD,stdout,PetscOptionsObject->prefix,next->option,next->text,next->man,next->flist,(char*)next->data,(char*)next->data);CHKERRQ(ierr);
3463cc1e11dSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3473cc1e11dSBarry Smith       if (str[0]) {
348e55864a3SBarry Smith         PetscOptionsObject->changedmethod = PETSC_TRUE;
3493cc1e11dSBarry Smith         next->set = PETSC_TRUE;
35064facd6cSBarry Smith         /* must use system malloc since SAWs may free this */
3515b02f95dSBarry Smith         ierr = PetscStrdup(str,(char**)&next->data);CHKERRQ(ierr);
3523cc1e11dSBarry Smith       }
3533cc1e11dSBarry Smith       break;
354b432afdaSMatthew Knepley     default:
355b432afdaSMatthew Knepley       break;
3566356e834SBarry Smith     }
3576356e834SBarry Smith     next = next->next;
3586356e834SBarry Smith   }
3596356e834SBarry Smith   PetscFunctionReturn(0);
3606356e834SBarry Smith }
3616356e834SBarry Smith 
362e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
363e04113cfSBarry Smith #include <petscviewersaws.h>
364d5649816SBarry Smith 
365d5649816SBarry Smith static int count = 0;
366d5649816SBarry Smith 
367e04113cfSBarry Smith PetscErrorCode PetscOptionsSAWsDestroy(void)
368d5649816SBarry Smith {
3692657e9d9SBarry Smith   PetscFunctionBegin;
370d5649816SBarry Smith   PetscFunctionReturn(0);
371d5649816SBarry Smith }
372d5649816SBarry Smith 
3739c1e0ce8SBarry Smith static const char *OptionsHeader = "<head>\n"
374a8d69d7bSBarry Smith                                    "<script type=\"text/javascript\" src=\"https://www.mcs.anl.gov/research/projects/saws/js/jquery-1.9.1.js\"></script>\n"
375a8d69d7bSBarry Smith                                    "<script type=\"text/javascript\" src=\"https://www.mcs.anl.gov/research/projects/saws/js/SAWs.js\"></script>\n"
376d1fc0251SBarry Smith                                    "<script type=\"text/javascript\" src=\"js/PETSc.js\"></script>\n"
37764bbc9efSBarry Smith                                    "<script>\n"
37864bbc9efSBarry Smith                                       "jQuery(document).ready(function() {\n"
37964bbc9efSBarry Smith                                          "PETSc.getAndDisplayDirectory(null,\"#variablesInfo\")\n"
38064bbc9efSBarry Smith                                       "})\n"
38164bbc9efSBarry Smith                                   "</script>\n"
38264bbc9efSBarry Smith                                   "</head>\n";
3831423471aSBarry Smith 
3841423471aSBarry Smith /*  Determines the size and style of the scroll region where PETSc options selectable from users are displayed */
3851423471aSBarry Smith static const char *OptionsBodyBottom = "<div id=\"variablesInfo\" style=\"background-color:lightblue;height:auto;max-height:500px;overflow:scroll;\"></div>\n<br>\n</body>";
38664bbc9efSBarry Smith 
387b3506946SBarry Smith /*
3887781c08eSBarry Smith     PetscOptionsSAWsInput - Presents all the PETSc Options processed by the program so the user may change them at runtime using the SAWs
389b3506946SBarry Smith 
390b3506946SBarry Smith     Bugs:
391b3506946SBarry Smith +    All processes must traverse through the exact same set of option queries do to the call to PetscScanString()
392b3506946SBarry Smith .    Internal strings have arbitrary length and string copies are not checked that they fit into string space
393b3506946SBarry Smith -    Only works for PetscInt == int, PetscReal == double etc
394b3506946SBarry Smith 
395b3506946SBarry Smith 
396b3506946SBarry Smith */
3974416b707SBarry Smith PetscErrorCode PetscOptionsSAWsInput(PetscOptionItems *PetscOptionsObject)
398b3506946SBarry Smith {
399b3506946SBarry Smith   PetscErrorCode  ierr;
4004416b707SBarry Smith   PetscOptionItem next     = PetscOptionsObject->next;
401d5649816SBarry Smith   static int      mancount = 0;
402b3506946SBarry Smith   char            options[16];
4037aab2a10SBarry Smith   PetscBool       changedmethod = PETSC_FALSE;
404a23eb982SSurtai Han   PetscBool       stopasking    = PETSC_FALSE;
40588a9752cSBarry Smith   char            manname[16],textname[16];
4062657e9d9SBarry Smith   char            dir[1024];
407b3506946SBarry Smith 
4082a409bb0SBarry Smith   PetscFunctionBegin;
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 
4127325c4a4SBarry Smith   PetscOptionsObject->pprefix = PetscOptionsObject->prefix; /* SAWs 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;
50912655325SBarry 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 
5174416b707SBarry Smith PetscErrorCode PetscOptionsEnd_Private(PetscOptionItems *PetscOptionsObject)
51853acd3b1SBarry Smith {
51953acd3b1SBarry Smith   PetscErrorCode  ierr;
5204416b707SBarry Smith   PetscOptionItem last;
5216356e834SBarry Smith   char            option[256],value[1024],tmp[32];
522330cf3c9SBarry Smith   size_t          j;
52353acd3b1SBarry Smith 
52453acd3b1SBarry Smith   PetscFunctionBegin;
52583355fc5SBarry Smith   if (PetscOptionsObject->next) {
52683355fc5SBarry Smith     if (!PetscOptionsObject->count) {
527a264d7a6SBarry Smith #if defined(PETSC_HAVE_SAWS)
528f7b25cf6SBarry Smith       ierr = PetscOptionsSAWsInput(PetscOptionsObject);CHKERRQ(ierr);
529b3506946SBarry Smith #else
530e55864a3SBarry Smith       ierr = PetscOptionsGetFromTextInput(PetscOptionsObject);CHKERRQ(ierr);
531b3506946SBarry Smith #endif
532aee2cecaSBarry Smith     }
533aee2cecaSBarry Smith   }
5346356e834SBarry Smith 
535e55864a3SBarry Smith   ierr = PetscFree(PetscOptionsObject->title);CHKERRQ(ierr);
5366356e834SBarry Smith 
537e26ddf31SBarry Smith   /* reset counter to -2; this updates the screen with the new options for the selected method */
538e55864a3SBarry Smith   if (PetscOptionsObject->changedmethod) PetscOptionsObject->count = -2;
5397a72a596SBarry Smith   /* reset alreadyprinted flag */
540e55864a3SBarry Smith   PetscOptionsObject->alreadyprinted = PETSC_FALSE;
541e55864a3SBarry Smith   if (PetscOptionsObject->object) PetscOptionsObject->object->optionsprinted = PETSC_TRUE;
542e55864a3SBarry Smith   PetscOptionsObject->object = NULL;
54353acd3b1SBarry Smith 
544e55864a3SBarry Smith   while (PetscOptionsObject->next) {
545e55864a3SBarry Smith     if (PetscOptionsObject->next->set) {
546e55864a3SBarry Smith       if (PetscOptionsObject->prefix) {
54753acd3b1SBarry Smith         ierr = PetscStrcpy(option,"-");CHKERRQ(ierr);
548e55864a3SBarry Smith         ierr = PetscStrcat(option,PetscOptionsObject->prefix);CHKERRQ(ierr);
549e55864a3SBarry Smith         ierr = PetscStrcat(option,PetscOptionsObject->next->option+1);CHKERRQ(ierr);
5506356e834SBarry Smith       } else {
551e55864a3SBarry Smith         ierr = PetscStrcpy(option,PetscOptionsObject->next->option);CHKERRQ(ierr);
5526356e834SBarry Smith       }
5536356e834SBarry Smith 
554e55864a3SBarry Smith       switch (PetscOptionsObject->next->type) {
5556356e834SBarry Smith       case OPTION_HEAD:
5566356e834SBarry Smith         break;
5576356e834SBarry Smith       case OPTION_INT_ARRAY:
558e55864a3SBarry Smith         sprintf(value,"%d",(int)((PetscInt*)PetscOptionsObject->next->data)[0]);
559e55864a3SBarry Smith         for (j=1; j<PetscOptionsObject->next->arraylength; j++) {
560e55864a3SBarry Smith           sprintf(tmp,"%d",(int)((PetscInt*)PetscOptionsObject->next->data)[j]);
5616356e834SBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5626356e834SBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5636356e834SBarry Smith         }
5646356e834SBarry Smith         break;
5656356e834SBarry Smith       case OPTION_INT:
566e55864a3SBarry Smith         sprintf(value,"%d",(int) *(PetscInt*)PetscOptionsObject->next->data);
5676356e834SBarry Smith         break;
5686356e834SBarry Smith       case OPTION_REAL:
569e55864a3SBarry Smith         sprintf(value,"%g",(double) *(PetscReal*)PetscOptionsObject->next->data);
5706356e834SBarry Smith         break;
5716356e834SBarry Smith       case OPTION_REAL_ARRAY:
572e55864a3SBarry Smith         sprintf(value,"%g",(double)((PetscReal*)PetscOptionsObject->next->data)[0]);
573e55864a3SBarry Smith         for (j=1; j<PetscOptionsObject->next->arraylength; j++) {
574e55864a3SBarry Smith           sprintf(tmp,"%g",(double)((PetscReal*)PetscOptionsObject->next->data)[j]);
5756356e834SBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5766356e834SBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5776356e834SBarry Smith         }
5786356e834SBarry Smith         break;
579050cccc3SHong Zhang       case OPTION_SCALAR_ARRAY:
58095f3a755SJose E. Roman         sprintf(value,"%g+%gi",(double)PetscRealPart(((PetscScalar*)PetscOptionsObject->next->data)[0]),(double)PetscImaginaryPart(((PetscScalar*)PetscOptionsObject->next->data)[0]));
581050cccc3SHong Zhang         for (j=1; j<PetscOptionsObject->next->arraylength; j++) {
58295f3a755SJose E. Roman           sprintf(tmp,"%g+%gi",(double)PetscRealPart(((PetscScalar*)PetscOptionsObject->next->data)[j]),(double)PetscImaginaryPart(((PetscScalar*)PetscOptionsObject->next->data)[j]));
583050cccc3SHong Zhang           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
584050cccc3SHong Zhang           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
585050cccc3SHong Zhang         }
586050cccc3SHong Zhang         break;
5877781c08eSBarry Smith       case OPTION_BOOL:
588e55864a3SBarry Smith         sprintf(value,"%d",*(int*)PetscOptionsObject->next->data);
5896356e834SBarry Smith         break;
5907781c08eSBarry Smith       case OPTION_BOOL_ARRAY:
591e55864a3SBarry Smith         sprintf(value,"%d",(int)((PetscBool*)PetscOptionsObject->next->data)[0]);
592e55864a3SBarry Smith         for (j=1; j<PetscOptionsObject->next->arraylength; j++) {
593e55864a3SBarry Smith           sprintf(tmp,"%d",(int)((PetscBool*)PetscOptionsObject->next->data)[j]);
5941ae3d29cSBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5951ae3d29cSBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5961ae3d29cSBarry Smith         }
5971ae3d29cSBarry Smith         break;
598a264d7a6SBarry Smith       case OPTION_FLIST:
5996991f827SBarry Smith         ierr = PetscStrcpy(value,(char*)PetscOptionsObject->next->data);CHKERRQ(ierr);
6006991f827SBarry Smith         break;
6011ae3d29cSBarry Smith       case OPTION_ELIST:
602e55864a3SBarry Smith         ierr = PetscStrcpy(value,(char*)PetscOptionsObject->next->data);CHKERRQ(ierr);
6036356e834SBarry Smith         break;
6041ae3d29cSBarry Smith       case OPTION_STRING:
605e55864a3SBarry Smith         ierr = PetscStrcpy(value,(char*)PetscOptionsObject->next->data);CHKERRQ(ierr);
606d51da6bfSBarry Smith         break;
6071ae3d29cSBarry Smith       case OPTION_STRING_ARRAY:
608e55864a3SBarry Smith         sprintf(value,"%s",((char**)PetscOptionsObject->next->data)[0]);
609e55864a3SBarry Smith         for (j=1; j<PetscOptionsObject->next->arraylength; j++) {
610e55864a3SBarry Smith           sprintf(tmp,"%s",((char**)PetscOptionsObject->next->data)[j]);
6111ae3d29cSBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
6121ae3d29cSBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
6131ae3d29cSBarry Smith         }
6146356e834SBarry Smith         break;
6156356e834SBarry Smith       }
616c5929fdfSBarry Smith       ierr = PetscOptionsSetValue(PetscOptionsObject->options,option,value);CHKERRQ(ierr);
6176356e834SBarry Smith     }
6186991f827SBarry Smith     if (PetscOptionsObject->next->type == OPTION_ELIST) {
6196991f827SBarry Smith       ierr = PetscStrNArrayDestroy(PetscOptionsObject->next->nlist,(char ***)&PetscOptionsObject->next->list);CHKERRQ(ierr);
6206991f827SBarry Smith     }
621e55864a3SBarry Smith     ierr   = PetscFree(PetscOptionsObject->next->text);CHKERRQ(ierr);
622e55864a3SBarry Smith     ierr   = PetscFree(PetscOptionsObject->next->option);CHKERRQ(ierr);
623e55864a3SBarry Smith     ierr   = PetscFree(PetscOptionsObject->next->man);CHKERRQ(ierr);
624e55864a3SBarry Smith     ierr   = PetscFree(PetscOptionsObject->next->edata);CHKERRQ(ierr);
625c979a496SBarry Smith 
62683355fc5SBarry Smith     if ((PetscOptionsObject->next->type == OPTION_STRING) || (PetscOptionsObject->next->type == OPTION_FLIST) || (PetscOptionsObject->next->type == OPTION_ELIST)){
62783355fc5SBarry Smith       free(PetscOptionsObject->next->data);
628c979a496SBarry Smith     } else {
62983355fc5SBarry Smith       ierr   = PetscFree(PetscOptionsObject->next->data);CHKERRQ(ierr);
630c979a496SBarry Smith     }
6317781c08eSBarry Smith 
63283355fc5SBarry Smith     last                    = PetscOptionsObject->next;
63383355fc5SBarry Smith     PetscOptionsObject->next = PetscOptionsObject->next->next;
6346356e834SBarry Smith     ierr                    = PetscFree(last);CHKERRQ(ierr);
6356356e834SBarry Smith   }
636f59f755dSBarry Smith   ierr = PetscFree(PetscOptionsObject->prefix);CHKERRQ(ierr);
637e55864a3SBarry Smith   PetscOptionsObject->next = 0;
63853acd3b1SBarry Smith   PetscFunctionReturn(0);
63953acd3b1SBarry Smith }
64053acd3b1SBarry Smith 
64153acd3b1SBarry Smith /*@C
64253acd3b1SBarry Smith    PetscOptionsEnum - Gets the enum value for a particular option in the database.
64353acd3b1SBarry Smith 
6443f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
64553acd3b1SBarry Smith 
64653acd3b1SBarry Smith    Input Parameters:
64753acd3b1SBarry Smith +  opt - option name
64853acd3b1SBarry Smith .  text - short string that describes the option
64953acd3b1SBarry Smith .  man - manual page with additional information on option
65053acd3b1SBarry Smith .  list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null
6510fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
6520fdccdaeSBarry Smith $                 PetscOptionsEnum(..., obj->value,&object->value,...) or
6530fdccdaeSBarry Smith $                 value = defaultvalue
6540fdccdaeSBarry Smith $                 PetscOptionsEnum(..., value,&value,&flg);
6550fdccdaeSBarry Smith $                 if (flg) {
65653acd3b1SBarry Smith 
65753acd3b1SBarry Smith    Output Parameter:
65853acd3b1SBarry Smith +  value - the  value to return
659b32e0204SMatthew G Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
66053acd3b1SBarry Smith 
66153acd3b1SBarry Smith    Level: beginner
66253acd3b1SBarry Smith 
66395452b02SPatrick Sanan    Notes:
66495452b02SPatrick Sanan     Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
66553acd3b1SBarry Smith 
66653acd3b1SBarry Smith           list is usually something like PCASMTypes or some other predefined list of enum names
66753acd3b1SBarry Smith 
6682efd9cb1SBarry Smith           If the user does not supply the option at all value is NOT changed. Thus
6692efd9cb1SBarry Smith           you should ALWAYS initialize value if you access it without first checking if the set flag is true.
6702efd9cb1SBarry Smith 
671989712b9SBarry Smith           The default/currentvalue passed into this routine does not get transferred to the output value variable automatically.
672989712b9SBarry Smith 
67353acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
674acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
675acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
67653acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
67753acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
678acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
679a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
68053acd3b1SBarry Smith @*/
6814416b707SBarry Smith PetscErrorCode  PetscOptionsEnum_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],const char *const *list,PetscEnum currentvalue,PetscEnum *value,PetscBool  *set)
68253acd3b1SBarry Smith {
68353acd3b1SBarry Smith   PetscErrorCode ierr;
68453acd3b1SBarry Smith   PetscInt       ntext = 0;
685aa5bb8c0SSatish Balay   PetscInt       tval;
686ace3abfcSBarry Smith   PetscBool      tflg;
68753acd3b1SBarry Smith 
68853acd3b1SBarry Smith   PetscFunctionBegin;
68953acd3b1SBarry Smith   while (list[ntext++]) {
690e32f2f54SBarry Smith     if (ntext > 50) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument appears to be wrong or have more than 50 entries");
69153acd3b1SBarry Smith   }
692e32f2f54SBarry Smith   if (ntext < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument must have at least two entries: typename and type prefix");
69353acd3b1SBarry Smith   ntext -= 3;
694e55864a3SBarry Smith   ierr   = PetscOptionsEList_Private(PetscOptionsObject,opt,text,man,list,ntext,list[currentvalue],&tval,&tflg);CHKERRQ(ierr);
695aa5bb8c0SSatish Balay   /* with PETSC_USE_64BIT_INDICES sizeof(PetscInt) != sizeof(PetscEnum) */
696aa5bb8c0SSatish Balay   if (tflg) *value = (PetscEnum)tval;
697aa5bb8c0SSatish Balay   if (set)  *set   = tflg;
69853acd3b1SBarry Smith   PetscFunctionReturn(0);
69953acd3b1SBarry Smith }
70053acd3b1SBarry Smith 
701d3e47460SLisandro Dalcin /*@C
702d3e47460SLisandro Dalcin    PetscOptionsEnumArray - Gets an array of enum values for a particular
703d3e47460SLisandro Dalcin    option in the database.
704d3e47460SLisandro Dalcin 
705d3e47460SLisandro Dalcin    Logically Collective on the communicator passed in PetscOptionsBegin()
706d3e47460SLisandro Dalcin 
707d3e47460SLisandro Dalcin    Input Parameters:
708d3e47460SLisandro Dalcin +  opt - the option one is seeking
709d3e47460SLisandro Dalcin .  text - short string describing option
710d3e47460SLisandro Dalcin .  man - manual page for option
71122d58ec6SMatthew G. Knepley .  list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null
712d3e47460SLisandro Dalcin -  n - maximum number of values
713d3e47460SLisandro Dalcin 
714d3e47460SLisandro Dalcin    Output Parameter:
715d3e47460SLisandro Dalcin +  value - location to copy values
716d3e47460SLisandro Dalcin .  n - actual number of values found
717d3e47460SLisandro Dalcin -  set - PETSC_TRUE if found, else PETSC_FALSE
718d3e47460SLisandro Dalcin 
719d3e47460SLisandro Dalcin    Level: beginner
720d3e47460SLisandro Dalcin 
721d3e47460SLisandro Dalcin    Notes:
722d3e47460SLisandro Dalcin    The array must be passed as a comma separated list.
723d3e47460SLisandro Dalcin 
724d3e47460SLisandro Dalcin    There must be no intervening spaces between the values.
725d3e47460SLisandro Dalcin 
726d3e47460SLisandro Dalcin    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
727d3e47460SLisandro Dalcin 
728d3e47460SLisandro Dalcin .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
729d3e47460SLisandro Dalcin           PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
730d3e47460SLisandro Dalcin           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
731d3e47460SLisandro Dalcin           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
732d3e47460SLisandro Dalcin           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
733d3e47460SLisandro Dalcin           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsRealArray()
734d3e47460SLisandro Dalcin @*/
7354416b707SBarry Smith PetscErrorCode  PetscOptionsEnumArray_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],const char *const *list,PetscEnum value[],PetscInt *n,PetscBool  *set)
736d3e47460SLisandro Dalcin {
737d3e47460SLisandro Dalcin   PetscInt        i,nlist = 0;
7384416b707SBarry Smith   PetscOptionItem amsopt;
739d3e47460SLisandro Dalcin   PetscErrorCode  ierr;
740d3e47460SLisandro Dalcin 
741d3e47460SLisandro Dalcin   PetscFunctionBegin;
742d3e47460SLisandro Dalcin   while (list[nlist++]) if (nlist > 50) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument appears to be wrong or have more than 50 entries");
743d3e47460SLisandro Dalcin   if (nlist < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument must have at least two entries: typename and type prefix");
744d3e47460SLisandro Dalcin   nlist -= 3; /* drop enum name, prefix, and null termination */
745d3e47460SLisandro Dalcin   if (0 && !PetscOptionsObject->count) { /* XXX Requires additional support */
746d3e47460SLisandro Dalcin     PetscEnum *vals;
7474416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_INT_ARRAY/*XXX OPTION_ENUM_ARRAY*/,&amsopt);CHKERRQ(ierr);
748d3e47460SLisandro Dalcin     ierr = PetscStrNArrayallocpy(nlist,list,(char***)&amsopt->list);CHKERRQ(ierr);
749d3e47460SLisandro Dalcin     amsopt->nlist = nlist;
750d3e47460SLisandro Dalcin     ierr = PetscMalloc1(*n,(PetscEnum**)&amsopt->data);CHKERRQ(ierr);
751d3e47460SLisandro Dalcin     amsopt->arraylength = *n;
752d3e47460SLisandro Dalcin     vals = (PetscEnum*)amsopt->data;
753d3e47460SLisandro Dalcin     for (i=0; i<*n; i++) vals[i] = value[i];
754d3e47460SLisandro Dalcin   }
755c5929fdfSBarry Smith   ierr = PetscOptionsGetEnumArray(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,list,value,n,set);CHKERRQ(ierr);
756d3e47460SLisandro Dalcin   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
757d3e47460SLisandro Dalcin     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%s",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,list[value[0]]);CHKERRQ(ierr);
758d3e47460SLisandro Dalcin     for (i=1; i<*n; i++) {ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%s",list[value[i]]);CHKERRQ(ierr);}
759d3e47460SLisandro Dalcin     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (choose from)",text);CHKERRQ(ierr);
760d3e47460SLisandro Dalcin     for (i=0; i<nlist; i++) {ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," %s",list[i]);CHKERRQ(ierr);}
761d3e47460SLisandro Dalcin     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," (%s)\n",ManSection(man));CHKERRQ(ierr);
762d3e47460SLisandro Dalcin   }
763d3e47460SLisandro Dalcin   PetscFunctionReturn(0);
764d3e47460SLisandro Dalcin }
765d3e47460SLisandro Dalcin 
766*5a856986SBarry Smith /*@C
767*5a856986SBarry Smith    PetscOptionsBoundedInt - Gets an integer value greater than or equal a given bound for a particular option in the database.
768*5a856986SBarry Smith 
769*5a856986SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
770*5a856986SBarry Smith 
771*5a856986SBarry Smith    Input Parameters:
772*5a856986SBarry Smith +  opt - option name
773*5a856986SBarry Smith .  text - short string that describes the option
774*5a856986SBarry Smith .  man - manual page with additional information on option
775*5a856986SBarry Smith .  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
776*5a856986SBarry Smith $                 PetscOptionsInt(..., obj->value,&object->value,...) or
777*5a856986SBarry Smith $                 value = defaultvalue
778*5a856986SBarry Smith $                 PetscOptionsInt(..., value,&value,&flg);
779*5a856986SBarry Smith $                 if (flg) {
780*5a856986SBarry Smith -  bound - the requested value should be large than this bound
781*5a856986SBarry Smith 
782*5a856986SBarry Smith    Output Parameter:
783*5a856986SBarry Smith +  value - the integer value to return
784*5a856986SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
785*5a856986SBarry Smith 
786*5a856986SBarry Smith    Notes:
787*5a856986SBarry Smith     If the user does not supply the option at all value is NOT changed. Thus
788*5a856986SBarry Smith     you should ALWAYS initialize value if you access it without first checking if the set flag is true.
789*5a856986SBarry Smith 
790*5a856986SBarry Smith     The default/currentvalue passed into this routine does not get transferred to the output value variable automatically.
791*5a856986SBarry Smith 
792*5a856986SBarry Smith     Errors if the supplied integer does not satisfy the bound
793*5a856986SBarry Smith 
794*5a856986SBarry Smith     Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
795*5a856986SBarry Smith 
796*5a856986SBarry Smith    Level: beginner
797*5a856986SBarry Smith 
798*5a856986SBarry Smith .seealso: PetscOptionsInt(), PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
799*5a856986SBarry Smith           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), PetscOptionsRangeInt()
800*5a856986SBarry Smith           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
801*5a856986SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
802*5a856986SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
803*5a856986SBarry Smith           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
804*5a856986SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
805*5a856986SBarry Smith @*/
806*5a856986SBarry Smith 
807*5a856986SBarry Smith /*@C
808*5a856986SBarry Smith    PetscOptionsRangeInt - Gets an integer value within a range of values for a particular option in the database.
809*5a856986SBarry Smith 
810*5a856986SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
811*5a856986SBarry Smith 
812*5a856986SBarry Smith    Input Parameters:
813*5a856986SBarry Smith +  opt - option name
814*5a856986SBarry Smith .  text - short string that describes the option
815*5a856986SBarry Smith .  man - manual page with additional information on option
816*5a856986SBarry Smith .  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
817*5a856986SBarry Smith $                 PetscOptionsInt(..., obj->value,&object->value,...) or
818*5a856986SBarry Smith $                 value = defaultvalue
819*5a856986SBarry Smith $                 PetscOptionsInt(..., value,&value,&flg);
820*5a856986SBarry Smith $                 if (flg) {
821*5a856986SBarry Smith .  lb - the lower bound
822*5a856986SBarry Smith -  ub - the upper bound
823*5a856986SBarry Smith 
824*5a856986SBarry Smith    Output Parameter:
825*5a856986SBarry Smith +  value - the integer value to return
826*5a856986SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
827*5a856986SBarry Smith 
828*5a856986SBarry Smith    Notes:
829*5a856986SBarry Smith     If the user does not supply the option at all value is NOT changed. Thus
830*5a856986SBarry Smith     you should ALWAYS initialize value if you access it without first checking if the set flag is true.
831*5a856986SBarry Smith 
832*5a856986SBarry Smith     The default/currentvalue passed into this routine does not get transferred to the output value variable automatically.
833*5a856986SBarry Smith 
834*5a856986SBarry Smith     Errors if the supplied integer does not satisfy the range
835*5a856986SBarry Smith 
836*5a856986SBarry Smith     Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
837*5a856986SBarry Smith 
838*5a856986SBarry Smith    Level: beginner
839*5a856986SBarry Smith 
840*5a856986SBarry Smith .seealso: PetscOptionsInt(), PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
841*5a856986SBarry Smith           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), PetscOptionsBoundedInt()
842*5a856986SBarry Smith           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
843*5a856986SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
844*5a856986SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
845*5a856986SBarry Smith           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
846*5a856986SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
847*5a856986SBarry Smith @*/
848*5a856986SBarry Smith 
84953acd3b1SBarry Smith /*@C
85053acd3b1SBarry Smith    PetscOptionsInt - Gets the integer value for a particular option in the database.
85153acd3b1SBarry Smith 
8523f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
85353acd3b1SBarry Smith 
85453acd3b1SBarry Smith    Input Parameters:
85553acd3b1SBarry Smith +  opt - option name
85653acd3b1SBarry Smith .  text - short string that describes the option
85753acd3b1SBarry Smith .  man - manual page with additional information on option
8580fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
8590fdccdaeSBarry Smith $                 PetscOptionsInt(..., obj->value,&object->value,...) or
8600fdccdaeSBarry Smith $                 value = defaultvalue
8610fdccdaeSBarry Smith $                 PetscOptionsInt(..., value,&value,&flg);
8620fdccdaeSBarry Smith $                 if (flg) {
86353acd3b1SBarry Smith 
86453acd3b1SBarry Smith    Output Parameter:
86553acd3b1SBarry Smith +  value - the integer value to return
86653acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
86753acd3b1SBarry Smith 
86895452b02SPatrick Sanan    Notes:
86995452b02SPatrick Sanan     If the user does not supply the option at all value is NOT changed. Thus
8702efd9cb1SBarry Smith     you should ALWAYS initialize value if you access it without first checking if the set flag is true.
8712efd9cb1SBarry Smith 
872989712b9SBarry Smith     The default/currentvalue passed into this routine does not get transferred to the output value variable automatically.
873989712b9SBarry Smith 
87495452b02SPatrick Sanan     Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
87553acd3b1SBarry Smith 
876*5a856986SBarry Smith    Level: beginner
877*5a856986SBarry Smith 
878*5a856986SBarry Smith .seealso: PetscOptionsBoundedInt(), PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
879*5a856986SBarry Smith           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), PetscOptionsRangeInt()
880acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
88153acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
88253acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
883acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
884a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
88553acd3b1SBarry Smith @*/
886*5a856986SBarry Smith PetscErrorCode  PetscOptionsInt_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscInt currentvalue,PetscInt *value,PetscBool  *set,PetscInt lb,PetscInt ub)
88753acd3b1SBarry Smith {
88853acd3b1SBarry Smith   PetscErrorCode  ierr;
8894416b707SBarry Smith   PetscOptionItem amsopt;
89012655325SBarry Smith   PetscBool       wasset;
89153acd3b1SBarry Smith 
89253acd3b1SBarry Smith   PetscFunctionBegin;
893*5a856986SBarry Smith   if (currentvalue < lb) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Current value %D less than allowed bound %D",currentvalue,lb);
894*5a856986SBarry Smith   if (currentvalue > ub ) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Newly set value %D greater than allowed bound %D",currentvalue,ub);
895e55864a3SBarry Smith      if (!PetscOptionsObject->count) {
8964416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_INT,&amsopt);CHKERRQ(ierr);
8976356e834SBarry Smith     ierr = PetscMalloc(sizeof(PetscInt),&amsopt->data);CHKERRQ(ierr);
89812655325SBarry Smith     *(PetscInt*)amsopt->data = currentvalue;
8993e211508SBarry Smith 
900c5929fdfSBarry Smith     ierr = PetscOptionsGetInt(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,&currentvalue,&wasset);CHKERRQ(ierr);
9013e211508SBarry Smith     if (wasset) {
90212655325SBarry Smith       *(PetscInt*)amsopt->data = currentvalue;
9033e211508SBarry Smith     }
904af6d86caSBarry Smith   }
90544ef3d73SBarry Smith   ierr = PetscOptionsGetInt(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,&wasset);CHKERRQ(ierr);
906*5a856986SBarry Smith   if (wasset && *value < lb ) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Newly set value %D less than allowed bound %D",*value,lb);
907*5a856986SBarry Smith   if (wasset && *value > ub ) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Newly set value %D greater than allowed bound %D",*value,ub);
90844ef3d73SBarry Smith   if (set) *set = wasset;
909e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
9100e4c290aSBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <now %D : formerly %D>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,wasset && value ? *value : currentvalue,currentvalue,text,ManSection(man));CHKERRQ(ierr);
91153acd3b1SBarry Smith   }
91253acd3b1SBarry Smith   PetscFunctionReturn(0);
91353acd3b1SBarry Smith }
91453acd3b1SBarry Smith 
91553acd3b1SBarry Smith /*@C
91653acd3b1SBarry Smith    PetscOptionsString - Gets the string value for a particular option in the database.
91753acd3b1SBarry Smith 
9183f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
91953acd3b1SBarry Smith 
92053acd3b1SBarry Smith    Input Parameters:
92153acd3b1SBarry Smith +  opt - option name
92253acd3b1SBarry Smith .  text - short string that describes the option
92353acd3b1SBarry Smith .  man - manual page with additional information on option
9240fdccdaeSBarry Smith .  currentvalue - the current value; caller is responsible for setting this value correctly. This is not used to set value
925bcbf2dc5SJed Brown -  len - length of the result string including null terminator
92653acd3b1SBarry Smith 
92753acd3b1SBarry Smith    Output Parameter:
92853acd3b1SBarry Smith +  value - the value to return
92953acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
93053acd3b1SBarry Smith 
93153acd3b1SBarry Smith    Level: beginner
93253acd3b1SBarry Smith 
93395452b02SPatrick Sanan    Notes:
93495452b02SPatrick Sanan     Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
93553acd3b1SBarry Smith 
9367fccdfe4SBarry 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).
9377fccdfe4SBarry Smith 
9382efd9cb1SBarry Smith           If the user does not supply the option at all value is NOT changed. Thus
9392efd9cb1SBarry Smith           you should ALWAYS initialize value if you access it without first checking if the set flag is true.
9402efd9cb1SBarry Smith 
941989712b9SBarry Smith           The default/currentvalue passed into this routine does not get transferred to the output value variable automatically.
942989712b9SBarry Smith 
9432efd9cb1SBarry Smith 
94489a13869SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
945acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
946acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsReal(), PetscOptionsBool(),
94753acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
94853acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
949acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
950a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
95153acd3b1SBarry Smith @*/
9524416b707SBarry Smith PetscErrorCode  PetscOptionsString_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],const char currentvalue[],char value[],size_t len,PetscBool  *set)
95353acd3b1SBarry Smith {
95453acd3b1SBarry Smith   PetscErrorCode  ierr;
9554416b707SBarry Smith   PetscOptionItem amsopt;
95644ef3d73SBarry Smith   PetscBool       lset;
95753acd3b1SBarry Smith 
95853acd3b1SBarry Smith   PetscFunctionBegin;
9591a1499c8SBarry Smith   if (!PetscOptionsObject->count) {
9604416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
96164facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
9620fdccdaeSBarry Smith     ierr = PetscStrdup(currentvalue ? currentvalue : "",(char**)&amsopt->data);CHKERRQ(ierr);
963af6d86caSBarry Smith   }
96444ef3d73SBarry Smith   ierr = PetscOptionsGetString(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,len,&lset);CHKERRQ(ierr);
96544ef3d73SBarry Smith   if (set) *set = lset;
966e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
9670e4c290aSBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <now %s : formerly %s>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,lset && value ? value : currentvalue,currentvalue,text,ManSection(man));CHKERRQ(ierr);
96853acd3b1SBarry Smith   }
96953acd3b1SBarry Smith   PetscFunctionReturn(0);
97053acd3b1SBarry Smith }
97153acd3b1SBarry Smith 
97253acd3b1SBarry Smith /*@C
97353acd3b1SBarry Smith    PetscOptionsReal - Gets the PetscReal value for a particular option in the database.
97453acd3b1SBarry Smith 
9753f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
97653acd3b1SBarry Smith 
97753acd3b1SBarry Smith    Input Parameters:
97853acd3b1SBarry Smith +  opt - option name
97953acd3b1SBarry Smith .  text - short string that describes the option
98053acd3b1SBarry Smith .  man - manual page with additional information on option
9810fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
9820fdccdaeSBarry Smith $                 PetscOptionsReal(..., obj->value,&object->value,...) or
9830fdccdaeSBarry Smith $                 value = defaultvalue
9840fdccdaeSBarry Smith $                 PetscOptionsReal(..., value,&value,&flg);
9850fdccdaeSBarry Smith $                 if (flg) {
98653acd3b1SBarry Smith 
98753acd3b1SBarry Smith    Output Parameter:
98853acd3b1SBarry Smith +  value - the value to return
98953acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
99053acd3b1SBarry Smith 
99195452b02SPatrick Sanan    Notes:
99295452b02SPatrick Sanan     If the user does not supply the option at all value is NOT changed. Thus
9932efd9cb1SBarry Smith     you should ALWAYS initialize value if you access it without first checking if the set flag is true.
9942efd9cb1SBarry Smith 
995989712b9SBarry Smith     The default/currentvalue passed into this routine does not get transferred to the output value variable automatically.
996989712b9SBarry Smith 
99795452b02SPatrick Sanan     Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
99853acd3b1SBarry Smith 
999*5a856986SBarry Smith    Level: beginner
1000*5a856986SBarry Smith 
100189a13869SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
1002acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
1003acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
100453acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
100553acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1006acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1007a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
100853acd3b1SBarry Smith @*/
10094416b707SBarry Smith PetscErrorCode  PetscOptionsReal_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscReal currentvalue,PetscReal *value,PetscBool  *set)
101053acd3b1SBarry Smith {
101153acd3b1SBarry Smith   PetscErrorCode  ierr;
10124416b707SBarry Smith   PetscOptionItem amsopt;
101344ef3d73SBarry Smith   PetscBool       lset;
101453acd3b1SBarry Smith 
101553acd3b1SBarry Smith   PetscFunctionBegin;
1016e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
10174416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_REAL,&amsopt);CHKERRQ(ierr);
1018538aa990SBarry Smith     ierr = PetscMalloc(sizeof(PetscReal),&amsopt->data);CHKERRQ(ierr);
1019a297a907SKarl Rupp 
10200fdccdaeSBarry Smith     *(PetscReal*)amsopt->data = currentvalue;
1021538aa990SBarry Smith   }
102244ef3d73SBarry Smith   ierr = PetscOptionsGetReal(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,&lset);CHKERRQ(ierr);
102344ef3d73SBarry Smith   if (set) *set = lset;
10241a1499c8SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
10250e4c290aSBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%g : %g>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,lset && value ? (double)*value : (double) currentvalue,(double)currentvalue,text,ManSection(man));CHKERRQ(ierr);
102653acd3b1SBarry Smith   }
102753acd3b1SBarry Smith   PetscFunctionReturn(0);
102853acd3b1SBarry Smith }
102953acd3b1SBarry Smith 
103053acd3b1SBarry Smith /*@C
103153acd3b1SBarry Smith    PetscOptionsScalar - Gets the scalar value for a particular option in the database.
103253acd3b1SBarry Smith 
10333f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
103453acd3b1SBarry Smith 
103553acd3b1SBarry Smith    Input Parameters:
103653acd3b1SBarry Smith +  opt - option name
103753acd3b1SBarry Smith .  text - short string that describes the option
103853acd3b1SBarry Smith .  man - manual page with additional information on option
10390fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
10400fdccdaeSBarry Smith $                 PetscOptionsScalar(..., obj->value,&object->value,...) or
10410fdccdaeSBarry Smith $                 value = defaultvalue
10420fdccdaeSBarry Smith $                 PetscOptionsScalar(..., value,&value,&flg);
10430fdccdaeSBarry Smith $                 if (flg) {
10440fdccdaeSBarry Smith 
104553acd3b1SBarry Smith 
104653acd3b1SBarry Smith    Output Parameter:
104753acd3b1SBarry Smith +  value - the value to return
104853acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
104953acd3b1SBarry Smith 
105095452b02SPatrick Sanan    Notes:
105195452b02SPatrick Sanan     If the user does not supply the option at all value is NOT changed. Thus
10522efd9cb1SBarry Smith     you should ALWAYS initialize value if you access it without first checking if the set flag is true.
10532efd9cb1SBarry Smith 
1054989712b9SBarry Smith     The default/currentvalue passed into this routine does not get transferred to the output value variable automatically.
1055989712b9SBarry Smith 
105695452b02SPatrick Sanan     Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
105753acd3b1SBarry Smith 
1058*5a856986SBarry Smith    Level: beginner
1059*5a856986SBarry Smith 
106089a13869SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
1061acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
1062acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
106353acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
106453acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1065acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1066a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
106753acd3b1SBarry Smith @*/
10684416b707SBarry Smith PetscErrorCode  PetscOptionsScalar_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscScalar currentvalue,PetscScalar *value,PetscBool  *set)
106953acd3b1SBarry Smith {
107053acd3b1SBarry Smith   PetscErrorCode ierr;
107153acd3b1SBarry Smith 
107253acd3b1SBarry Smith   PetscFunctionBegin;
107353acd3b1SBarry Smith #if !defined(PETSC_USE_COMPLEX)
10740fdccdaeSBarry Smith   ierr = PetscOptionsReal(opt,text,man,currentvalue,value,set);CHKERRQ(ierr);
107553acd3b1SBarry Smith #else
1076c5929fdfSBarry Smith   ierr = PetscOptionsGetScalar(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,set);CHKERRQ(ierr);
107753acd3b1SBarry Smith #endif
107853acd3b1SBarry Smith   PetscFunctionReturn(0);
107953acd3b1SBarry Smith }
108053acd3b1SBarry Smith 
108153acd3b1SBarry Smith /*@C
108290d69ab7SBarry 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
108390d69ab7SBarry Smith                       its value is set to false.
108453acd3b1SBarry Smith 
10853f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
108653acd3b1SBarry Smith 
108753acd3b1SBarry Smith    Input Parameters:
108853acd3b1SBarry Smith +  opt - option name
108953acd3b1SBarry Smith .  text - short string that describes the option
109053acd3b1SBarry Smith -  man - manual page with additional information on option
109153acd3b1SBarry Smith 
109253acd3b1SBarry Smith    Output Parameter:
109353acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
109453acd3b1SBarry Smith 
109553acd3b1SBarry Smith    Level: beginner
109653acd3b1SBarry Smith 
109795452b02SPatrick Sanan    Notes:
109895452b02SPatrick Sanan     Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
109953acd3b1SBarry Smith 
110089a13869SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
1101acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
1102acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
110353acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
110453acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1105acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1106a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
110753acd3b1SBarry Smith @*/
11084416b707SBarry Smith PetscErrorCode  PetscOptionsName_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool  *flg)
110953acd3b1SBarry Smith {
111053acd3b1SBarry Smith   PetscErrorCode  ierr;
11114416b707SBarry Smith   PetscOptionItem amsopt;
111253acd3b1SBarry Smith 
111353acd3b1SBarry Smith   PetscFunctionBegin;
1114e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
11154416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1116ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1117a297a907SKarl Rupp 
1118ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
11191ae3d29cSBarry Smith   }
1120c5929fdfSBarry Smith   ierr = PetscOptionsHasName(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,flg);CHKERRQ(ierr);
1121e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1122e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
112353acd3b1SBarry Smith   }
112453acd3b1SBarry Smith   PetscFunctionReturn(0);
112553acd3b1SBarry Smith }
112653acd3b1SBarry Smith 
112753acd3b1SBarry Smith /*@C
1128a264d7a6SBarry Smith      PetscOptionsFList - Puts a list of option values that a single one may be selected from
112953acd3b1SBarry Smith 
11303f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
113153acd3b1SBarry Smith 
113253acd3b1SBarry Smith    Input Parameters:
113353acd3b1SBarry Smith +  opt - option name
113453acd3b1SBarry Smith .  text - short string that describes the option
113553acd3b1SBarry Smith .  man - manual page with additional information on option
113653acd3b1SBarry Smith .  list - the possible choices
11370fdccdaeSBarry Smith .  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with
11380fdccdaeSBarry Smith $                 PetscOptionsFlist(..., obj->value,value,len,&flg);
11390fdccdaeSBarry Smith $                 if (flg) {
11403cc1e11dSBarry Smith -  len - the length of the character array value
114153acd3b1SBarry Smith 
114253acd3b1SBarry Smith    Output Parameter:
114353acd3b1SBarry Smith +  value - the value to return
114453acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
114553acd3b1SBarry Smith 
114653acd3b1SBarry Smith    Level: intermediate
114753acd3b1SBarry Smith 
114895452b02SPatrick Sanan    Notes:
114995452b02SPatrick Sanan     Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
115053acd3b1SBarry Smith 
11512efd9cb1SBarry Smith           If the user does not supply the option at all value is NOT changed. Thus
11522efd9cb1SBarry Smith           you should ALWAYS initialize value if you access it without first checking if the set flag is true.
11532efd9cb1SBarry Smith 
1154989712b9SBarry Smith           The default/currentvalue passed into this routine does not get transferred to the output value variable automatically.
1155989712b9SBarry Smith 
115653acd3b1SBarry Smith    See PetscOptionsEList() for when the choices are given in a string array
115753acd3b1SBarry Smith 
115853acd3b1SBarry Smith    To get a listing of all currently specified options,
115988c29154SBarry Smith     see PetscOptionsView() or PetscOptionsGetAll()
116053acd3b1SBarry Smith 
1161eabe10d7SBarry Smith    Developer Note: This cannot check for invalid selection because of things like MATAIJ that are not included in the list
1162eabe10d7SBarry Smith 
116389a13869SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1164acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
116553acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
116653acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1167acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1168a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsEnum()
116953acd3b1SBarry Smith @*/
11704416b707SBarry Smith PetscErrorCode  PetscOptionsFList_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char ltext[],const char man[],PetscFunctionList list,const char currentvalue[],char value[],size_t len,PetscBool  *set)
117153acd3b1SBarry Smith {
117253acd3b1SBarry Smith   PetscErrorCode  ierr;
11734416b707SBarry Smith   PetscOptionItem amsopt;
117444ef3d73SBarry Smith   PetscBool       lset;
117553acd3b1SBarry Smith 
117653acd3b1SBarry Smith   PetscFunctionBegin;
11771a1499c8SBarry Smith   if (!PetscOptionsObject->count) {
11784416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,ltext,man,OPTION_FLIST,&amsopt);CHKERRQ(ierr);
117964facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
11800fdccdaeSBarry Smith     ierr = PetscStrdup(currentvalue ? currentvalue : "",(char**)&amsopt->data);CHKERRQ(ierr);
11813cc1e11dSBarry Smith     amsopt->flist = list;
11823cc1e11dSBarry Smith   }
118344ef3d73SBarry Smith   ierr = PetscOptionsGetString(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,len,&lset);CHKERRQ(ierr);
118444ef3d73SBarry Smith   if (set) *set = lset;
11851a1499c8SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
118644ef3d73SBarry Smith     ierr = PetscFunctionListPrintTypes(PetscOptionsObject->comm,stdout,PetscOptionsObject->prefix,opt,ltext,man,list,currentvalue,lset && value ? value : currentvalue);CHKERRQ(ierr);CHKERRQ(ierr);
118753acd3b1SBarry Smith   }
118853acd3b1SBarry Smith   PetscFunctionReturn(0);
118953acd3b1SBarry Smith }
119053acd3b1SBarry Smith 
119153acd3b1SBarry Smith /*@C
119253acd3b1SBarry Smith      PetscOptionsEList - Puts a list of option values that a single one may be selected from
119353acd3b1SBarry Smith 
11943f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
119553acd3b1SBarry Smith 
119653acd3b1SBarry Smith    Input Parameters:
119753acd3b1SBarry Smith +  opt - option name
119853acd3b1SBarry Smith .  ltext - short string that describes the option
119953acd3b1SBarry Smith .  man - manual page with additional information on option
1200a264d7a6SBarry Smith .  list - the possible choices (one of these must be selected, anything else is invalid)
120153acd3b1SBarry Smith .  ntext - number of choices
12020fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with
12030fdccdaeSBarry Smith $                 PetscOptionsElist(..., obj->value,&value,&flg);
12040fdccdaeSBarry Smith $                 if (flg) {
12050fdccdaeSBarry Smith 
120653acd3b1SBarry Smith 
120753acd3b1SBarry Smith    Output Parameter:
120853acd3b1SBarry Smith +  value - the index of the value to return
120953acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
121053acd3b1SBarry Smith 
121153acd3b1SBarry Smith    Level: intermediate
121253acd3b1SBarry Smith 
121395452b02SPatrick Sanan    Notes:
121495452b02SPatrick Sanan     Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
121553acd3b1SBarry Smith 
12162efd9cb1SBarry Smith          If the user does not supply the option at all value is NOT changed. Thus
12172efd9cb1SBarry Smith           you should ALWAYS initialize value if you access it without first checking if the set flag is true.
12182efd9cb1SBarry Smith 
1219a264d7a6SBarry Smith    See PetscOptionsFList() for when the choices are given in a PetscFunctionList()
122053acd3b1SBarry Smith 
122189a13869SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1222acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
122353acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
122453acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1225acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1226a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEnum()
122753acd3b1SBarry Smith @*/
12284416b707SBarry Smith PetscErrorCode  PetscOptionsEList_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char ltext[],const char man[],const char *const *list,PetscInt ntext,const char currentvalue[],PetscInt *value,PetscBool  *set)
122953acd3b1SBarry Smith {
123053acd3b1SBarry Smith   PetscErrorCode  ierr;
123153acd3b1SBarry Smith   PetscInt        i;
12324416b707SBarry Smith   PetscOptionItem amsopt;
123344ef3d73SBarry Smith   PetscBool       lset;
123453acd3b1SBarry Smith 
123553acd3b1SBarry Smith   PetscFunctionBegin;
12361a1499c8SBarry Smith   if (!PetscOptionsObject->count) {
12374416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,ltext,man,OPTION_ELIST,&amsopt);CHKERRQ(ierr);
123864facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
12390fdccdaeSBarry Smith     ierr = PetscStrdup(currentvalue ? currentvalue : "",(char**)&amsopt->data);CHKERRQ(ierr);
12406991f827SBarry Smith     ierr = PetscStrNArrayallocpy(ntext,list,(char***)&amsopt->list);CHKERRQ(ierr);
12411ae3d29cSBarry Smith     amsopt->nlist = ntext;
12421ae3d29cSBarry Smith   }
124344ef3d73SBarry Smith   ierr = PetscOptionsGetEList(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,list,ntext,value,&lset);CHKERRQ(ierr);
124444ef3d73SBarry Smith   if (set) *set = lset;
12451a1499c8SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
12460e4c290aSBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <now %s : formerly %s> %s (choose one of)",PetscOptionsObject->prefix?PetscOptionsObject->prefix:"",opt+1,lset && value ? list[*value] : currentvalue,currentvalue,ltext);CHKERRQ(ierr);
124753acd3b1SBarry Smith     for (i=0; i<ntext; i++) {
1248e55864a3SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," %s",list[i]);CHKERRQ(ierr);
124953acd3b1SBarry Smith     }
1250e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," (%s)\n",ManSection(man));CHKERRQ(ierr);
125153acd3b1SBarry Smith   }
125253acd3b1SBarry Smith   PetscFunctionReturn(0);
125353acd3b1SBarry Smith }
125453acd3b1SBarry Smith 
125553acd3b1SBarry Smith /*@C
1256acfcf0e5SJed Brown      PetscOptionsBoolGroupBegin - First in a series of logical queries on the options database for
1257d5649816SBarry Smith        which at most a single value can be true.
125853acd3b1SBarry Smith 
12593f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
126053acd3b1SBarry Smith 
126153acd3b1SBarry Smith    Input Parameters:
126253acd3b1SBarry Smith +  opt - option name
126353acd3b1SBarry Smith .  text - short string that describes the option
126453acd3b1SBarry Smith -  man - manual page with additional information on option
126553acd3b1SBarry Smith 
126653acd3b1SBarry Smith    Output Parameter:
126753acd3b1SBarry Smith .  flg - whether that option was set or not
126853acd3b1SBarry Smith 
126953acd3b1SBarry Smith    Level: intermediate
127053acd3b1SBarry Smith 
127195452b02SPatrick Sanan    Notes:
127295452b02SPatrick Sanan     Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
127353acd3b1SBarry Smith 
1274acfcf0e5SJed Brown    Must be followed by 0 or more PetscOptionsBoolGroup()s and PetscOptionsBoolGroupEnd()
127553acd3b1SBarry Smith 
127689a13869SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1277acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
127853acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
127953acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1280acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1281a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
128253acd3b1SBarry Smith @*/
12834416b707SBarry Smith PetscErrorCode  PetscOptionsBoolGroupBegin_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool  *flg)
128453acd3b1SBarry Smith {
128553acd3b1SBarry Smith   PetscErrorCode  ierr;
12864416b707SBarry Smith   PetscOptionItem amsopt;
128753acd3b1SBarry Smith 
128853acd3b1SBarry Smith   PetscFunctionBegin;
1289e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
12904416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1291ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1292a297a907SKarl Rupp 
1293ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
12941ae3d29cSBarry Smith   }
129568b16fdaSBarry Smith   *flg = PETSC_FALSE;
1296c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,flg,NULL);CHKERRQ(ierr);
1297e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1298e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  Pick at most one of -------------\n");CHKERRQ(ierr);
1299e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"    -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
130053acd3b1SBarry Smith   }
130153acd3b1SBarry Smith   PetscFunctionReturn(0);
130253acd3b1SBarry Smith }
130353acd3b1SBarry Smith 
130453acd3b1SBarry Smith /*@C
1305acfcf0e5SJed Brown      PetscOptionsBoolGroup - One in a series of logical queries on the options database for
1306d5649816SBarry Smith        which at most a single value can be true.
130753acd3b1SBarry Smith 
13083f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
130953acd3b1SBarry Smith 
131053acd3b1SBarry Smith    Input Parameters:
131153acd3b1SBarry Smith +  opt - option name
131253acd3b1SBarry Smith .  text - short string that describes the option
131353acd3b1SBarry Smith -  man - manual page with additional information on option
131453acd3b1SBarry Smith 
131553acd3b1SBarry Smith    Output Parameter:
131653acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
131753acd3b1SBarry Smith 
131853acd3b1SBarry Smith    Level: intermediate
131953acd3b1SBarry Smith 
132095452b02SPatrick Sanan    Notes:
132195452b02SPatrick Sanan     Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
132253acd3b1SBarry Smith 
1323acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin() and preceded a PetscOptionsBoolGroupEnd()
132453acd3b1SBarry Smith 
132589a13869SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1326acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
132753acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
132853acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1329acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1330a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
133153acd3b1SBarry Smith @*/
13324416b707SBarry Smith PetscErrorCode  PetscOptionsBoolGroup_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool  *flg)
133353acd3b1SBarry Smith {
133453acd3b1SBarry Smith   PetscErrorCode  ierr;
13354416b707SBarry Smith   PetscOptionItem amsopt;
133653acd3b1SBarry Smith 
133753acd3b1SBarry Smith   PetscFunctionBegin;
1338e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
13394416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1340ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1341a297a907SKarl Rupp 
1342ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
13431ae3d29cSBarry Smith   }
134417326d04SJed Brown   *flg = PETSC_FALSE;
1345c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,flg,NULL);CHKERRQ(ierr);
1346e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1347e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"    -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
134853acd3b1SBarry Smith   }
134953acd3b1SBarry Smith   PetscFunctionReturn(0);
135053acd3b1SBarry Smith }
135153acd3b1SBarry Smith 
135253acd3b1SBarry Smith /*@C
1353acfcf0e5SJed Brown      PetscOptionsBoolGroupEnd - Last in a series of logical queries on the options database for
1354d5649816SBarry Smith        which at most a single value can be true.
135553acd3b1SBarry Smith 
13563f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
135753acd3b1SBarry Smith 
135853acd3b1SBarry Smith    Input Parameters:
135953acd3b1SBarry Smith +  opt - option name
136053acd3b1SBarry Smith .  text - short string that describes the option
136153acd3b1SBarry Smith -  man - manual page with additional information on option
136253acd3b1SBarry Smith 
136353acd3b1SBarry Smith    Output Parameter:
136453acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
136553acd3b1SBarry Smith 
136653acd3b1SBarry Smith    Level: intermediate
136753acd3b1SBarry Smith 
136895452b02SPatrick Sanan    Notes:
136995452b02SPatrick Sanan     Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
137053acd3b1SBarry Smith 
1371acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin()
137253acd3b1SBarry Smith 
137389a13869SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1374acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
137553acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
137653acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1377acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1378a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
137953acd3b1SBarry Smith @*/
13804416b707SBarry Smith PetscErrorCode  PetscOptionsBoolGroupEnd_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool  *flg)
138153acd3b1SBarry Smith {
138253acd3b1SBarry Smith   PetscErrorCode  ierr;
13834416b707SBarry Smith   PetscOptionItem amsopt;
138453acd3b1SBarry Smith 
138553acd3b1SBarry Smith   PetscFunctionBegin;
1386e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
13874416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1388ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1389a297a907SKarl Rupp 
1390ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
13911ae3d29cSBarry Smith   }
139217326d04SJed Brown   *flg = PETSC_FALSE;
1393c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,flg,NULL);CHKERRQ(ierr);
1394e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1395e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"    -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
139653acd3b1SBarry Smith   }
139753acd3b1SBarry Smith   PetscFunctionReturn(0);
139853acd3b1SBarry Smith }
139953acd3b1SBarry Smith 
140053acd3b1SBarry Smith /*@C
1401acfcf0e5SJed Brown    PetscOptionsBool - Determines if a particular option is in the database with a true or false
140253acd3b1SBarry Smith 
14033f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
140453acd3b1SBarry Smith 
140553acd3b1SBarry Smith    Input Parameters:
140653acd3b1SBarry Smith +  opt - option name
140753acd3b1SBarry Smith .  text - short string that describes the option
1408868c398cSBarry Smith .  man - manual page with additional information on option
140994ae4db5SBarry Smith -  currentvalue - the current value
141053acd3b1SBarry Smith 
141153acd3b1SBarry Smith    Output Parameter:
141253acd3b1SBarry Smith .  flg - PETSC_TRUE or PETSC_FALSE
141353acd3b1SBarry Smith .  set - PETSC_TRUE if found, else PETSC_FALSE
141453acd3b1SBarry Smith 
14152efd9cb1SBarry Smith    Notes:
14162efd9cb1SBarry Smith        TRUE, true, YES, yes, nostring, and 1 all translate to PETSC_TRUE
14172efd9cb1SBarry Smith        FALSE, false, NO, no, and 0 all translate to PETSC_FALSE
14182efd9cb1SBarry Smith 
14192efd9cb1SBarry Smith       If the option is given, but no value is provided, then flg and set are both given the value PETSC_TRUE. That is -requested_bool
14202efd9cb1SBarry Smith      is equivalent to -requested_bool true
14212efd9cb1SBarry Smith 
14222efd9cb1SBarry Smith        If the user does not supply the option at all flg is NOT changed. Thus
14232efd9cb1SBarry Smith      you should ALWAYS initialize the flg if you access it without first checking if the set flag is true.
14242efd9cb1SBarry Smith 
142595452b02SPatrick Sanan     Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
142653acd3b1SBarry Smith 
1427*5a856986SBarry Smith    Level: beginner
1428*5a856986SBarry Smith 
142989a13869SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
1430acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
1431acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
143253acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
143353acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1434acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1435a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
143653acd3b1SBarry Smith @*/
14374416b707SBarry Smith PetscErrorCode  PetscOptionsBool_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool currentvalue,PetscBool  *flg,PetscBool  *set)
143853acd3b1SBarry Smith {
143953acd3b1SBarry Smith   PetscErrorCode  ierr;
1440ace3abfcSBarry Smith   PetscBool       iset;
14414416b707SBarry Smith   PetscOptionItem amsopt;
144253acd3b1SBarry Smith 
144353acd3b1SBarry Smith   PetscFunctionBegin;
1444e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
14454416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1446ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1447a297a907SKarl Rupp 
144894ae4db5SBarry Smith     *(PetscBool*)amsopt->data = currentvalue;
1449af6d86caSBarry Smith   }
1450c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,flg,&iset);CHKERRQ(ierr);
145153acd3b1SBarry Smith   if (set) *set = iset;
14521a1499c8SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
145344ef3d73SBarry Smith     const char *v = PetscBools[currentvalue], *vn = PetscBools[iset && flg ? *flg : currentvalue];
145444ef3d73SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s: <%s : %s> %s (%s)\n",PetscOptionsObject->prefix?PetscOptionsObject->prefix:"",opt+1,v,vn,text,ManSection(man));CHKERRQ(ierr);
145553acd3b1SBarry Smith   }
145653acd3b1SBarry Smith   PetscFunctionReturn(0);
145753acd3b1SBarry Smith }
145853acd3b1SBarry Smith 
145953acd3b1SBarry Smith /*@C
146053acd3b1SBarry Smith    PetscOptionsRealArray - Gets an array of double values for a particular
146153acd3b1SBarry Smith    option in the database. The values must be separated with commas with
146253acd3b1SBarry Smith    no intervening spaces.
146353acd3b1SBarry Smith 
14643f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
146553acd3b1SBarry Smith 
146653acd3b1SBarry Smith    Input Parameters:
146753acd3b1SBarry Smith +  opt - the option one is seeking
146853acd3b1SBarry Smith .  text - short string describing option
146953acd3b1SBarry Smith .  man - manual page for option
147053acd3b1SBarry Smith -  nmax - maximum number of values
147153acd3b1SBarry Smith 
147253acd3b1SBarry Smith    Output Parameter:
147353acd3b1SBarry Smith +  value - location to copy values
147453acd3b1SBarry Smith .  nmax - actual number of values found
147553acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
147653acd3b1SBarry Smith 
147753acd3b1SBarry Smith    Level: beginner
147853acd3b1SBarry Smith 
147953acd3b1SBarry Smith    Notes:
148053acd3b1SBarry Smith    The user should pass in an array of doubles
148153acd3b1SBarry Smith 
148253acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
148353acd3b1SBarry Smith 
148489a13869SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1485acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
148653acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
148753acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1488acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1489a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
149053acd3b1SBarry Smith @*/
14914416b707SBarry Smith PetscErrorCode PetscOptionsRealArray_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscReal value[],PetscInt *n,PetscBool  *set)
149253acd3b1SBarry Smith {
149353acd3b1SBarry Smith   PetscErrorCode  ierr;
149453acd3b1SBarry Smith   PetscInt        i;
14954416b707SBarry Smith   PetscOptionItem amsopt;
149653acd3b1SBarry Smith 
149753acd3b1SBarry Smith   PetscFunctionBegin;
1498e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
1499e26ddf31SBarry Smith     PetscReal *vals;
1500e26ddf31SBarry Smith 
15014416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_REAL_ARRAY,&amsopt);CHKERRQ(ierr);
1502e55864a3SBarry Smith     ierr = PetscMalloc((*n)*sizeof(PetscReal),&amsopt->data);CHKERRQ(ierr);
1503e26ddf31SBarry Smith     vals = (PetscReal*)amsopt->data;
1504e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1505e26ddf31SBarry Smith     amsopt->arraylength = *n;
1506e26ddf31SBarry Smith   }
1507c5929fdfSBarry Smith   ierr = PetscOptionsGetRealArray(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,n,set);CHKERRQ(ierr);
1508e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1509a519f713SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%g",PetscOptionsObject->prefix?PetscOptionsObject->prefix:"",opt+1,(double)value[0]);CHKERRQ(ierr);
151053acd3b1SBarry Smith     for (i=1; i<*n; i++) {
1511e55864a3SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%g",(double)value[i]);CHKERRQ(ierr);
151253acd3b1SBarry Smith     }
1513e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
151453acd3b1SBarry Smith   }
151553acd3b1SBarry Smith   PetscFunctionReturn(0);
151653acd3b1SBarry Smith }
151753acd3b1SBarry Smith 
1518050cccc3SHong Zhang /*@C
1519050cccc3SHong Zhang    PetscOptionsScalarArray - Gets an array of Scalar values for a particular
1520050cccc3SHong Zhang    option in the database. The values must be separated with commas with
1521050cccc3SHong Zhang    no intervening spaces.
1522050cccc3SHong Zhang 
1523050cccc3SHong Zhang    Logically Collective on the communicator passed in PetscOptionsBegin()
1524050cccc3SHong Zhang 
1525050cccc3SHong Zhang    Input Parameters:
1526050cccc3SHong Zhang +  opt - the option one is seeking
1527050cccc3SHong Zhang .  text - short string describing option
1528050cccc3SHong Zhang .  man - manual page for option
1529050cccc3SHong Zhang -  nmax - maximum number of values
1530050cccc3SHong Zhang 
1531050cccc3SHong Zhang    Output Parameter:
1532050cccc3SHong Zhang +  value - location to copy values
1533050cccc3SHong Zhang .  nmax - actual number of values found
1534050cccc3SHong Zhang -  set - PETSC_TRUE if found, else PETSC_FALSE
1535050cccc3SHong Zhang 
1536050cccc3SHong Zhang    Level: beginner
1537050cccc3SHong Zhang 
1538050cccc3SHong Zhang    Notes:
1539050cccc3SHong Zhang    The user should pass in an array of doubles
1540050cccc3SHong Zhang 
1541050cccc3SHong Zhang    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
1542050cccc3SHong Zhang 
154389a13869SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1544050cccc3SHong Zhang           PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
1545050cccc3SHong Zhang           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1546050cccc3SHong Zhang           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1547050cccc3SHong Zhang           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1548050cccc3SHong Zhang           PetscOptionsFList(), PetscOptionsEList()
1549050cccc3SHong Zhang @*/
15504416b707SBarry Smith PetscErrorCode PetscOptionsScalarArray_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscScalar value[],PetscInt *n,PetscBool  *set)
1551050cccc3SHong Zhang {
1552050cccc3SHong Zhang   PetscErrorCode  ierr;
1553050cccc3SHong Zhang   PetscInt        i;
15544416b707SBarry Smith   PetscOptionItem amsopt;
1555050cccc3SHong Zhang 
1556050cccc3SHong Zhang   PetscFunctionBegin;
1557050cccc3SHong Zhang   if (!PetscOptionsObject->count) {
1558050cccc3SHong Zhang     PetscScalar *vals;
1559050cccc3SHong Zhang 
15604416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_SCALAR_ARRAY,&amsopt);CHKERRQ(ierr);
1561050cccc3SHong Zhang     ierr = PetscMalloc((*n)*sizeof(PetscScalar),&amsopt->data);CHKERRQ(ierr);
1562050cccc3SHong Zhang     vals = (PetscScalar*)amsopt->data;
1563050cccc3SHong Zhang     for (i=0; i<*n; i++) vals[i] = value[i];
1564050cccc3SHong Zhang     amsopt->arraylength = *n;
1565050cccc3SHong Zhang   }
1566c5929fdfSBarry Smith   ierr = PetscOptionsGetScalarArray(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,n,set);CHKERRQ(ierr);
1567050cccc3SHong Zhang   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
156895f3a755SJose E. Roman     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%g+%gi",PetscOptionsObject->prefix?PetscOptionsObject->prefix:"",opt+1,(double)PetscRealPart(value[0]),(double)PetscImaginaryPart(value[0]));CHKERRQ(ierr);
1569050cccc3SHong Zhang     for (i=1; i<*n; i++) {
157095f3a755SJose E. Roman       ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%g+%gi",(double)PetscRealPart(value[i]),(double)PetscImaginaryPart(value[i]));CHKERRQ(ierr);
1571050cccc3SHong Zhang     }
1572050cccc3SHong Zhang     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
1573050cccc3SHong Zhang   }
1574050cccc3SHong Zhang   PetscFunctionReturn(0);
1575050cccc3SHong Zhang }
157653acd3b1SBarry Smith 
157753acd3b1SBarry Smith /*@C
157853acd3b1SBarry Smith    PetscOptionsIntArray - Gets an array of integers for a particular
1579b32a342fSShri Abhyankar    option in the database.
158053acd3b1SBarry Smith 
15813f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
158253acd3b1SBarry Smith 
158353acd3b1SBarry Smith    Input Parameters:
158453acd3b1SBarry Smith +  opt - the option one is seeking
158553acd3b1SBarry Smith .  text - short string describing option
158653acd3b1SBarry Smith .  man - manual page for option
1587f8a50e2bSBarry Smith -  n - maximum number of values
158853acd3b1SBarry Smith 
158953acd3b1SBarry Smith    Output Parameter:
159053acd3b1SBarry Smith +  value - location to copy values
1591f8a50e2bSBarry Smith .  n - actual number of values found
159253acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
159353acd3b1SBarry Smith 
159453acd3b1SBarry Smith    Level: beginner
159553acd3b1SBarry Smith 
159653acd3b1SBarry Smith    Notes:
1597b32a342fSShri Abhyankar    The array can be passed as
1598bebe2cf6SSatish Balay    a comma separated list:                                 0,1,2,3,4,5,6,7
15990fd488f5SShri Abhyankar    a range (start-end+1):                                  0-8
16000fd488f5SShri Abhyankar    a range with given increment (start-end+1:inc):         0-7:2
1601bebe2cf6SSatish Balay    a combination of values and ranges separated by commas: 0,1-8,8-15:2
1602b32a342fSShri Abhyankar 
1603b32a342fSShri Abhyankar    There must be no intervening spaces between the values.
160453acd3b1SBarry Smith 
160553acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
160653acd3b1SBarry Smith 
160789a13869SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1608acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
160953acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
161053acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1611acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1612a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsRealArray()
161353acd3b1SBarry Smith @*/
16144416b707SBarry Smith PetscErrorCode  PetscOptionsIntArray_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscInt value[],PetscInt *n,PetscBool  *set)
161553acd3b1SBarry Smith {
161653acd3b1SBarry Smith   PetscErrorCode ierr;
161753acd3b1SBarry Smith   PetscInt        i;
16184416b707SBarry Smith   PetscOptionItem amsopt;
161953acd3b1SBarry Smith 
162053acd3b1SBarry Smith   PetscFunctionBegin;
1621e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
1622e26ddf31SBarry Smith     PetscInt *vals;
1623e26ddf31SBarry Smith 
16244416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_INT_ARRAY,&amsopt);CHKERRQ(ierr);
1625854ce69bSBarry Smith     ierr = PetscMalloc1(*n,(PetscInt**)&amsopt->data);CHKERRQ(ierr);
1626e26ddf31SBarry Smith     vals = (PetscInt*)amsopt->data;
1627e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1628e26ddf31SBarry Smith     amsopt->arraylength = *n;
1629e26ddf31SBarry Smith   }
1630c5929fdfSBarry Smith   ierr = PetscOptionsGetIntArray(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,n,set);CHKERRQ(ierr);
1631e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1632e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%d",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,value[0]);CHKERRQ(ierr);
163353acd3b1SBarry Smith     for (i=1; i<*n; i++) {
1634e55864a3SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%d",value[i]);CHKERRQ(ierr);
163553acd3b1SBarry Smith     }
1636e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
163753acd3b1SBarry Smith   }
163853acd3b1SBarry Smith   PetscFunctionReturn(0);
163953acd3b1SBarry Smith }
164053acd3b1SBarry Smith 
164153acd3b1SBarry Smith /*@C
164253acd3b1SBarry Smith    PetscOptionsStringArray - Gets an array of string values for a particular
164353acd3b1SBarry Smith    option in the database. The values must be separated with commas with
164453acd3b1SBarry Smith    no intervening spaces.
164553acd3b1SBarry Smith 
16463f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
164753acd3b1SBarry Smith 
164853acd3b1SBarry Smith    Input Parameters:
164953acd3b1SBarry Smith +  opt - the option one is seeking
165053acd3b1SBarry Smith .  text - short string describing option
165153acd3b1SBarry Smith .  man - manual page for option
165253acd3b1SBarry Smith -  nmax - maximum number of strings
165353acd3b1SBarry Smith 
165453acd3b1SBarry Smith    Output Parameter:
165553acd3b1SBarry Smith +  value - location to copy strings
165653acd3b1SBarry Smith .  nmax - actual number of strings found
165753acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
165853acd3b1SBarry Smith 
165953acd3b1SBarry Smith    Level: beginner
166053acd3b1SBarry Smith 
166153acd3b1SBarry Smith    Notes:
166253acd3b1SBarry Smith    The user should pass in an array of pointers to char, to hold all the
166353acd3b1SBarry Smith    strings returned by this function.
166453acd3b1SBarry Smith 
166553acd3b1SBarry Smith    The user is responsible for deallocating the strings that are
166653acd3b1SBarry Smith    returned. The Fortran interface for this routine is not supported.
166753acd3b1SBarry Smith 
166853acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
166953acd3b1SBarry Smith 
167089a13869SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1671acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
167253acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
167353acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1674acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1675a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
167653acd3b1SBarry Smith @*/
16774416b707SBarry Smith PetscErrorCode  PetscOptionsStringArray_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],char *value[],PetscInt *nmax,PetscBool  *set)
167853acd3b1SBarry Smith {
167953acd3b1SBarry Smith   PetscErrorCode  ierr;
16804416b707SBarry Smith   PetscOptionItem amsopt;
168153acd3b1SBarry Smith 
168253acd3b1SBarry Smith   PetscFunctionBegin;
1683e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
16844416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_STRING_ARRAY,&amsopt);CHKERRQ(ierr);
1685854ce69bSBarry Smith     ierr = PetscMalloc1(*nmax,(char**)&amsopt->data);CHKERRQ(ierr);
1686a297a907SKarl Rupp 
16871ae3d29cSBarry Smith     amsopt->arraylength = *nmax;
16881ae3d29cSBarry Smith   }
1689c5929fdfSBarry Smith   ierr = PetscOptionsGetStringArray(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,nmax,set);CHKERRQ(ierr);
1690e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1691e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <string1,string2,...>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
169253acd3b1SBarry Smith   }
169353acd3b1SBarry Smith   PetscFunctionReturn(0);
169453acd3b1SBarry Smith }
169553acd3b1SBarry Smith 
1696e2446a98SMatthew Knepley /*@C
1697acfcf0e5SJed Brown    PetscOptionsBoolArray - Gets an array of logical values (true or false) for a particular
1698e2446a98SMatthew Knepley    option in the database. The values must be separated with commas with
1699e2446a98SMatthew Knepley    no intervening spaces.
1700e2446a98SMatthew Knepley 
17013f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
1702e2446a98SMatthew Knepley 
1703e2446a98SMatthew Knepley    Input Parameters:
1704e2446a98SMatthew Knepley +  opt - the option one is seeking
1705e2446a98SMatthew Knepley .  text - short string describing option
1706e2446a98SMatthew Knepley .  man - manual page for option
1707e2446a98SMatthew Knepley -  nmax - maximum number of values
1708e2446a98SMatthew Knepley 
1709e2446a98SMatthew Knepley    Output Parameter:
1710e2446a98SMatthew Knepley +  value - location to copy values
1711e2446a98SMatthew Knepley .  nmax - actual number of values found
1712e2446a98SMatthew Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
1713e2446a98SMatthew Knepley 
1714e2446a98SMatthew Knepley    Level: beginner
1715e2446a98SMatthew Knepley 
1716e2446a98SMatthew Knepley    Notes:
1717e2446a98SMatthew Knepley    The user should pass in an array of doubles
1718e2446a98SMatthew Knepley 
1719e2446a98SMatthew Knepley    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
1720e2446a98SMatthew Knepley 
172189a13869SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1722acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
1723e2446a98SMatthew Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1724e2446a98SMatthew Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1725acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1726a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
1727e2446a98SMatthew Knepley @*/
17284416b707SBarry Smith PetscErrorCode  PetscOptionsBoolArray_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool value[],PetscInt *n,PetscBool *set)
1729e2446a98SMatthew Knepley {
1730e2446a98SMatthew Knepley   PetscErrorCode   ierr;
1731e2446a98SMatthew Knepley   PetscInt         i;
17324416b707SBarry Smith   PetscOptionItem  amsopt;
1733e2446a98SMatthew Knepley 
1734e2446a98SMatthew Knepley   PetscFunctionBegin;
1735e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
1736ace3abfcSBarry Smith     PetscBool *vals;
17371ae3d29cSBarry Smith 
17384416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL_ARRAY,&amsopt);CHKERRQ(ierr);
17391a1499c8SBarry Smith     ierr = PetscMalloc1(*n,(PetscBool**)&amsopt->data);CHKERRQ(ierr);
1740ace3abfcSBarry Smith     vals = (PetscBool*)amsopt->data;
17411ae3d29cSBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
17421ae3d29cSBarry Smith     amsopt->arraylength = *n;
17431ae3d29cSBarry Smith   }
1744c5929fdfSBarry Smith   ierr = PetscOptionsGetBoolArray(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,n,set);CHKERRQ(ierr);
1745e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1746e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%d",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,value[0]);CHKERRQ(ierr);
1747e2446a98SMatthew Knepley     for (i=1; i<*n; i++) {
1748e55864a3SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%d",value[i]);CHKERRQ(ierr);
1749e2446a98SMatthew Knepley     }
1750e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
1751e2446a98SMatthew Knepley   }
1752e2446a98SMatthew Knepley   PetscFunctionReturn(0);
1753e2446a98SMatthew Knepley }
1754e2446a98SMatthew Knepley 
17558cc676e6SMatthew G Knepley /*@C
1756d1da0b69SBarry Smith    PetscOptionsViewer - Gets a viewer appropriate for the type indicated by the user
17578cc676e6SMatthew G Knepley 
17588cc676e6SMatthew G Knepley    Logically Collective on the communicator passed in PetscOptionsBegin()
17598cc676e6SMatthew G Knepley 
17608cc676e6SMatthew G Knepley    Input Parameters:
17618cc676e6SMatthew G Knepley +  opt - option name
17628cc676e6SMatthew G Knepley .  text - short string that describes the option
17638cc676e6SMatthew G Knepley -  man - manual page with additional information on option
17648cc676e6SMatthew G Knepley 
17658cc676e6SMatthew G Knepley    Output Parameter:
17668cc676e6SMatthew G Knepley +  viewer - the viewer
17678cc676e6SMatthew G Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
17688cc676e6SMatthew G Knepley 
17698cc676e6SMatthew G Knepley    Level: beginner
17708cc676e6SMatthew G Knepley 
177195452b02SPatrick Sanan    Notes:
177295452b02SPatrick Sanan     Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
17738cc676e6SMatthew G Knepley 
17745a7113b9SPatrick Sanan    See PetscOptionsGetViewer() for the format of the supplied viewer and its options
17758cc676e6SMatthew G Knepley 
177689a13869SBarry Smith .seealso: PetscOptionsGetViewer(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
17778cc676e6SMatthew G Knepley           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
17788cc676e6SMatthew G Knepley           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
17798cc676e6SMatthew G Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
17808cc676e6SMatthew G Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
17818cc676e6SMatthew G Knepley           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1782a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
17838cc676e6SMatthew G Knepley @*/
17844416b707SBarry Smith PetscErrorCode  PetscOptionsViewer_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscViewer *viewer,PetscViewerFormat *format,PetscBool  *set)
17858cc676e6SMatthew G Knepley {
17868cc676e6SMatthew G Knepley   PetscErrorCode  ierr;
17874416b707SBarry Smith   PetscOptionItem amsopt;
17888cc676e6SMatthew G Knepley 
17898cc676e6SMatthew G Knepley   PetscFunctionBegin;
17901a1499c8SBarry Smith   if (!PetscOptionsObject->count) {
17914416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
179264facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
17935b02f95dSBarry Smith     ierr = PetscStrdup("",(char**)&amsopt->data);CHKERRQ(ierr);
17948cc676e6SMatthew G Knepley   }
179516413a6aSBarry Smith   ierr = PetscOptionsGetViewer(PetscOptionsObject->comm,PetscOptionsObject->options,PetscOptionsObject->prefix,opt,viewer,format,set);CHKERRQ(ierr);
1796e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1797e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%s>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,"",text,ManSection(man));CHKERRQ(ierr);
17988cc676e6SMatthew G Knepley   }
17998cc676e6SMatthew G Knepley   PetscFunctionReturn(0);
18008cc676e6SMatthew G Knepley }
18018cc676e6SMatthew G Knepley 
180253acd3b1SBarry Smith 
180353acd3b1SBarry Smith /*@C
1804b52f573bSBarry Smith      PetscOptionsHead - Puts a heading before listing any more published options. Used, for example,
180553acd3b1SBarry Smith             in KSPSetFromOptions_GMRES().
180653acd3b1SBarry Smith 
18073f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
180853acd3b1SBarry Smith 
180953acd3b1SBarry Smith    Input Parameter:
181053acd3b1SBarry Smith .   head - the heading text
181153acd3b1SBarry Smith 
181253acd3b1SBarry Smith 
181353acd3b1SBarry Smith    Level: intermediate
181453acd3b1SBarry Smith 
181595452b02SPatrick Sanan    Notes:
181695452b02SPatrick Sanan     Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
181753acd3b1SBarry Smith 
1818b52f573bSBarry Smith           Can be followed by a call to PetscOptionsTail() in the same function.
181953acd3b1SBarry Smith 
182089a13869SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1821acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
182253acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
182353acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1824acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1825a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
182653acd3b1SBarry Smith @*/
18274416b707SBarry Smith PetscErrorCode  PetscOptionsHead(PetscOptionItems *PetscOptionsObject,const char head[])
182853acd3b1SBarry Smith {
182953acd3b1SBarry Smith   PetscErrorCode ierr;
183053acd3b1SBarry Smith 
183153acd3b1SBarry Smith   PetscFunctionBegin;
1832e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1833e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  %s\n",head);CHKERRQ(ierr);
183453acd3b1SBarry Smith   }
183553acd3b1SBarry Smith   PetscFunctionReturn(0);
183653acd3b1SBarry Smith }
183753acd3b1SBarry Smith 
183853acd3b1SBarry Smith 
183953acd3b1SBarry Smith 
184053acd3b1SBarry Smith 
184153acd3b1SBarry Smith 
184253acd3b1SBarry Smith 
1843