xref: /petsc/src/sys/objects/aoptions.c (revision 0fdccdaef8d180f33759be5ddfe586b8087048a7)
17d0a6c19SBarry Smith 
253acd3b1SBarry Smith /*
33fc1eb6aSBarry Smith    Implements the higher-level options database querying methods. These are self-documenting and can attach at runtime to
43fc1eb6aSBarry Smith    GUI code to display the options and get values from the users.
553acd3b1SBarry Smith 
653acd3b1SBarry Smith */
753acd3b1SBarry Smith 
8afcb2eb5SJed Brown #include <petsc-private/petscimpl.h>        /*I  "petscsys.h"   I*/
9665c2dedSJed Brown #include <petscviewer.h>
1053acd3b1SBarry Smith 
112aa6d131SJed Brown #define ManSection(str) ((str) ? (str) : "None")
122aa6d131SJed Brown 
1353acd3b1SBarry Smith /*
1453acd3b1SBarry Smith     Keep a linked list of options that have been posted and we are waiting for
153fc1eb6aSBarry Smith    user selection. See the manual page for PetscOptionsBegin()
1653acd3b1SBarry Smith 
1753acd3b1SBarry Smith     Eventually we'll attach this beast to a MPI_Comm
1853acd3b1SBarry Smith */
19f8d0b74dSMatthew Knepley PetscOptionsObjectType PetscOptionsObject;
2053acd3b1SBarry Smith PetscInt               PetscOptionsPublishCount = 0;
2153acd3b1SBarry Smith 
2253acd3b1SBarry Smith #undef __FUNCT__
2353acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsBegin_Private"
2453acd3b1SBarry Smith /*
2553acd3b1SBarry Smith     Handles setting up the data structure in a call to PetscOptionsBegin()
2653acd3b1SBarry Smith */
2753acd3b1SBarry Smith PetscErrorCode PetscOptionsBegin_Private(MPI_Comm comm,const char prefix[],const char title[],const char mansec[])
2853acd3b1SBarry Smith {
2953acd3b1SBarry Smith   PetscErrorCode ierr;
3053acd3b1SBarry Smith 
3153acd3b1SBarry Smith   PetscFunctionBegin;
3253acd3b1SBarry Smith   PetscOptionsObject.next          = 0;
3353acd3b1SBarry Smith   PetscOptionsObject.comm          = comm;
346356e834SBarry Smith   PetscOptionsObject.changedmethod = PETSC_FALSE;
35a297a907SKarl Rupp 
36c31cb41cSBarry Smith   ierr = PetscFree(PetscOptionsObject.prefix);CHKERRQ(ierr);
3753acd3b1SBarry Smith   ierr = PetscStrallocpy(prefix,&PetscOptionsObject.prefix);CHKERRQ(ierr);
38c31cb41cSBarry Smith   ierr = PetscFree(PetscOptionsObject.title);CHKERRQ(ierr);
3953acd3b1SBarry Smith   ierr = PetscStrallocpy(title,&PetscOptionsObject.title);CHKERRQ(ierr);
4053acd3b1SBarry Smith 
410298fd71SBarry Smith   ierr = PetscOptionsHasName(NULL,"-help",&PetscOptionsObject.printhelp);CHKERRQ(ierr);
4253acd3b1SBarry Smith   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1) {
4361b37b28SSatish Balay     if (!PetscOptionsObject.alreadyprinted) {
4453acd3b1SBarry Smith       ierr = (*PetscHelpPrintf)(comm,"%s -------------------------------------------------\n",title);CHKERRQ(ierr);
4553acd3b1SBarry Smith     }
4661b37b28SSatish Balay   }
4753acd3b1SBarry Smith   PetscFunctionReturn(0);
4853acd3b1SBarry Smith }
4953acd3b1SBarry Smith 
503194b578SJed Brown #undef __FUNCT__
513194b578SJed Brown #define __FUNCT__ "PetscObjectOptionsBegin_Private"
523194b578SJed Brown /*
533194b578SJed Brown     Handles setting up the data structure in a call to PetscObjectOptionsBegin()
543194b578SJed Brown */
553194b578SJed Brown PetscErrorCode PetscObjectOptionsBegin_Private(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);
633194b578SJed Brown   PetscOptionsObject.object         = obj;
643194b578SJed Brown   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   }
723194b578SJed Brown   ierr = PetscOptionsBegin_Private(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 */
7953acd3b1SBarry Smith #undef __FUNCT__
8053acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsCreate_Private"
81e3ed6ec8SBarry Smith static int PetscOptionsCreate_Private(const char opt[],const char text[],const char man[],PetscOptionType t,PetscOptions *amsopt)
8253acd3b1SBarry Smith {
8353acd3b1SBarry Smith   int          ierr;
8453acd3b1SBarry Smith   PetscOptions next;
853be6e4c3SJed Brown   PetscBool    valid;
8653acd3b1SBarry Smith 
8753acd3b1SBarry Smith   PetscFunctionBegin;
883be6e4c3SJed Brown   ierr = PetscOptionsValidKey(opt,&valid);CHKERRQ(ierr);
893be6e4c3SJed Brown   if (!valid) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"The option '%s' is not a valid key",opt);
903be6e4c3SJed Brown 
91b00a9115SJed Brown   ierr            = PetscNew(amsopt);CHKERRQ(ierr);
9253acd3b1SBarry Smith   (*amsopt)->next = 0;
9353acd3b1SBarry Smith   (*amsopt)->set  = PETSC_FALSE;
946356e834SBarry Smith   (*amsopt)->type = t;
9553acd3b1SBarry Smith   (*amsopt)->data = 0;
9661b37b28SSatish Balay 
9753acd3b1SBarry Smith   ierr = PetscStrallocpy(text,&(*amsopt)->text);CHKERRQ(ierr);
9853acd3b1SBarry Smith   ierr = PetscStrallocpy(opt,&(*amsopt)->option);CHKERRQ(ierr);
996356e834SBarry Smith   ierr = PetscStrallocpy(man,&(*amsopt)->man);CHKERRQ(ierr);
10053acd3b1SBarry Smith 
101a297a907SKarl Rupp   if (!PetscOptionsObject.next) PetscOptionsObject.next = *amsopt;
102a297a907SKarl Rupp   else {
10353acd3b1SBarry Smith     next = PetscOptionsObject.next;
10453acd3b1SBarry Smith     while (next->next) next = next->next;
10553acd3b1SBarry Smith     next->next = *amsopt;
10653acd3b1SBarry Smith   }
10753acd3b1SBarry Smith   PetscFunctionReturn(0);
10853acd3b1SBarry Smith }
10953acd3b1SBarry Smith 
11053acd3b1SBarry Smith #undef __FUNCT__
111aee2cecaSBarry Smith #define __FUNCT__ "PetscScanString"
112aee2cecaSBarry Smith /*
1133fc1eb6aSBarry Smith     PetscScanString -  Gets user input via stdin from process and broadcasts to all processes
1143fc1eb6aSBarry Smith 
1153fc1eb6aSBarry Smith     Collective on MPI_Comm
1163fc1eb6aSBarry Smith 
1173fc1eb6aSBarry Smith    Input Parameters:
1183fc1eb6aSBarry Smith +     commm - communicator for the broadcast, must be PETSC_COMM_WORLD
1193fc1eb6aSBarry Smith .     n - length of the string, must be the same on all processes
1203fc1eb6aSBarry Smith -     str - location to store input
121aee2cecaSBarry Smith 
122aee2cecaSBarry Smith     Bugs:
123aee2cecaSBarry Smith .   Assumes process 0 of the given communicator has access to stdin
124aee2cecaSBarry Smith 
125aee2cecaSBarry Smith */
1263fc1eb6aSBarry Smith static PetscErrorCode PetscScanString(MPI_Comm comm,size_t n,char str[])
127aee2cecaSBarry Smith {
128330cf3c9SBarry Smith   size_t         i;
129aee2cecaSBarry Smith   char           c;
1303fc1eb6aSBarry Smith   PetscMPIInt    rank,nm;
131aee2cecaSBarry Smith   PetscErrorCode ierr;
132aee2cecaSBarry Smith 
133aee2cecaSBarry Smith   PetscFunctionBegin;
134aee2cecaSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
135aee2cecaSBarry Smith   if (!rank) {
136aee2cecaSBarry Smith     c = (char) getchar();
137aee2cecaSBarry Smith     i = 0;
138aee2cecaSBarry Smith     while (c != '\n' && i < n-1) {
139aee2cecaSBarry Smith       str[i++] = c;
140aee2cecaSBarry Smith       c = (char)getchar();
141aee2cecaSBarry Smith     }
142aee2cecaSBarry Smith     str[i] = 0;
143aee2cecaSBarry Smith   }
1444dc2109aSBarry Smith   ierr = PetscMPIIntCast(n,&nm);CHKERRQ(ierr);
1453fc1eb6aSBarry Smith   ierr = MPI_Bcast(str,nm,MPI_CHAR,0,comm);CHKERRQ(ierr);
146aee2cecaSBarry Smith   PetscFunctionReturn(0);
147aee2cecaSBarry Smith }
148aee2cecaSBarry Smith 
149aca04275SShao-Ching Huang #undef __FUNCT__
150aca04275SShao-Ching Huang #define __FUNCT__ "PetscStrdup"
1515b02f95dSBarry Smith /*
1525b02f95dSBarry Smith     This is needed because certain strings may be freed by SAWs, hence we cannot use PetscStrallocpy()
1535b02f95dSBarry Smith */
1545b02f95dSBarry Smith static PetscErrorCode  PetscStrdup(const char s[],char *t[])
1555b02f95dSBarry Smith {
1565b02f95dSBarry Smith   PetscErrorCode ierr;
1575b02f95dSBarry Smith   size_t         len;
1585b02f95dSBarry Smith   char           *tmp = 0;
1595b02f95dSBarry Smith 
1605b02f95dSBarry Smith   PetscFunctionBegin;
1615b02f95dSBarry Smith   if (s) {
1625b02f95dSBarry Smith     ierr = PetscStrlen(s,&len);CHKERRQ(ierr);
163ee48884fSBarry Smith     tmp = (char*) malloc((len+1)*sizeof(char*));
1645b02f95dSBarry Smith     if (!tmp) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"No memory to duplicate string");
1655b02f95dSBarry Smith     ierr = PetscStrcpy(tmp,s);CHKERRQ(ierr);
1665b02f95dSBarry Smith   }
1675b02f95dSBarry Smith   *t = tmp;
1685b02f95dSBarry Smith   PetscFunctionReturn(0);
1695b02f95dSBarry Smith }
1705b02f95dSBarry Smith 
1715b02f95dSBarry Smith 
172aee2cecaSBarry Smith #undef __FUNCT__
173aee2cecaSBarry Smith #define __FUNCT__ "PetscOptionsGetFromTextInput"
174aee2cecaSBarry Smith /*
1753cc1e11dSBarry Smith     PetscOptionsGetFromTextInput - Presents all the PETSc Options processed by the program so the user may change them at runtime
176aee2cecaSBarry Smith 
177aee2cecaSBarry Smith     Notes: this isn't really practical, it is just to demonstrate the principle
178aee2cecaSBarry Smith 
1797781c08eSBarry Smith     A carriage return indicates no change from the default; but this like -ksp_monitor <stdout>  the default is actually not stdout the default
1807781c08eSBarry Smith     is to do nothing so to get it to use stdout you need to type stdout. This is kind of bug?
1817781c08eSBarry Smith 
182aee2cecaSBarry Smith     Bugs:
1837781c08eSBarry Smith +    All processes must traverse through the exact same set of option queries due to the call to PetscScanString()
1843cc1e11dSBarry Smith .    Internal strings have arbitrary length and string copies are not checked that they fit into string space
185aee2cecaSBarry Smith -    Only works for PetscInt == int, PetscReal == double etc
186aee2cecaSBarry Smith 
1873cc1e11dSBarry Smith     Developer Notes: Normally the GUI that presents the options the user and retrieves the values would be running in a different
1883cc1e11dSBarry Smith      address space and communicating with the PETSc program
1893cc1e11dSBarry Smith 
190aee2cecaSBarry Smith */
191aee2cecaSBarry Smith PetscErrorCode PetscOptionsGetFromTextInput()
1926356e834SBarry Smith {
1936356e834SBarry Smith   PetscErrorCode ierr;
1946356e834SBarry Smith   PetscOptions   next = PetscOptionsObject.next;
1956356e834SBarry Smith   char           str[512];
1967781c08eSBarry Smith   PetscBool      bid;
197a4404d99SBarry Smith   PetscReal      ir,*valr;
198330cf3c9SBarry Smith   PetscInt       *vald;
199330cf3c9SBarry Smith   size_t         i;
2006356e834SBarry Smith 
201e26ddf31SBarry Smith   ierr = (*PetscPrintf)(PETSC_COMM_WORLD,"%s -------------------------------------------------\n",PetscOptionsObject.title);CHKERRQ(ierr);
2026356e834SBarry Smith   while (next) {
2036356e834SBarry Smith     switch (next->type) {
2046356e834SBarry Smith     case OPTION_HEAD:
2056356e834SBarry Smith       break;
206e26ddf31SBarry Smith     case OPTION_INT_ARRAY:
207e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",next->option+1);CHKERRQ(ierr);
208e26ddf31SBarry Smith       vald = (PetscInt*) next->data;
209e26ddf31SBarry Smith       for (i=0; i<next->arraylength; i++) {
210e26ddf31SBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"%d",vald[i]);CHKERRQ(ierr);
211e26ddf31SBarry Smith         if (i < next->arraylength-1) {
212e26ddf31SBarry Smith           ierr = PetscPrintf(PETSC_COMM_WORLD,",");CHKERRQ(ierr);
213e26ddf31SBarry Smith         }
214e26ddf31SBarry Smith       }
215e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,">: %s (%s) ",next->text,next->man);CHKERRQ(ierr);
216e26ddf31SBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
217e26ddf31SBarry Smith       if (str[0]) {
218e26ddf31SBarry Smith         PetscToken token;
219e26ddf31SBarry Smith         PetscInt   n=0,nmax = next->arraylength,*dvalue = (PetscInt*)next->data,start,end;
220e26ddf31SBarry Smith         size_t     len;
221e26ddf31SBarry Smith         char       *value;
222ace3abfcSBarry Smith         PetscBool  foundrange;
223e26ddf31SBarry Smith 
224e26ddf31SBarry Smith         next->set = PETSC_TRUE;
225e26ddf31SBarry Smith         value     = str;
226e26ddf31SBarry Smith         ierr      = PetscTokenCreate(value,',',&token);CHKERRQ(ierr);
227e26ddf31SBarry Smith         ierr      = PetscTokenFind(token,&value);CHKERRQ(ierr);
228e26ddf31SBarry Smith         while (n < nmax) {
229e26ddf31SBarry Smith           if (!value) break;
230e26ddf31SBarry Smith 
231e26ddf31SBarry Smith           /* look for form  d-D where d and D are integers */
232e26ddf31SBarry Smith           foundrange = PETSC_FALSE;
233e26ddf31SBarry Smith           ierr       = PetscStrlen(value,&len);CHKERRQ(ierr);
234e26ddf31SBarry Smith           if (value[0] == '-') i=2;
235e26ddf31SBarry Smith           else i=1;
236330cf3c9SBarry Smith           for (;i<len; i++) {
237e26ddf31SBarry Smith             if (value[i] == '-') {
238e32f2f54SBarry Smith               if (i == len-1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Error in %D-th array entry %s\n",n,value);
239e26ddf31SBarry Smith               value[i] = 0;
240cfbddea1SSatish Balay               ierr     = PetscOptionsStringToInt(value,&start);CHKERRQ(ierr);
241cfbddea1SSatish Balay               ierr     = PetscOptionsStringToInt(value+i+1,&end);CHKERRQ(ierr);
242e32f2f54SBarry 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);
243e32f2f54SBarry 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);
244e26ddf31SBarry Smith               for (; start<end; start++) {
245e26ddf31SBarry Smith                 *dvalue = start; dvalue++;n++;
246e26ddf31SBarry Smith               }
247e26ddf31SBarry Smith               foundrange = PETSC_TRUE;
248e26ddf31SBarry Smith               break;
249e26ddf31SBarry Smith             }
250e26ddf31SBarry Smith           }
251e26ddf31SBarry Smith           if (!foundrange) {
252cfbddea1SSatish Balay             ierr = PetscOptionsStringToInt(value,dvalue);CHKERRQ(ierr);
253e26ddf31SBarry Smith             dvalue++;
254e26ddf31SBarry Smith             n++;
255e26ddf31SBarry Smith           }
256e26ddf31SBarry Smith           ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
257e26ddf31SBarry Smith         }
2588c74ee41SBarry Smith         ierr = PetscTokenDestroy(&token);CHKERRQ(ierr);
259e26ddf31SBarry Smith       }
260e26ddf31SBarry Smith       break;
261e26ddf31SBarry Smith     case OPTION_REAL_ARRAY:
262e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",next->option+1);CHKERRQ(ierr);
263e26ddf31SBarry Smith       valr = (PetscReal*) next->data;
264e26ddf31SBarry Smith       for (i=0; i<next->arraylength; i++) {
265e26ddf31SBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"%g",valr[i]);CHKERRQ(ierr);
266e26ddf31SBarry Smith         if (i < next->arraylength-1) {
267e26ddf31SBarry Smith           ierr = PetscPrintf(PETSC_COMM_WORLD,",");CHKERRQ(ierr);
268e26ddf31SBarry Smith         }
269e26ddf31SBarry Smith       }
270e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,">: %s (%s) ",next->text,next->man);CHKERRQ(ierr);
271e26ddf31SBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
272e26ddf31SBarry Smith       if (str[0]) {
273e26ddf31SBarry Smith         PetscToken token;
274e26ddf31SBarry Smith         PetscInt   n = 0,nmax = next->arraylength;
275e26ddf31SBarry Smith         PetscReal  *dvalue = (PetscReal*)next->data;
276e26ddf31SBarry Smith         char       *value;
277e26ddf31SBarry Smith 
278e26ddf31SBarry Smith         next->set = PETSC_TRUE;
279e26ddf31SBarry Smith         value     = str;
280e26ddf31SBarry Smith         ierr      = PetscTokenCreate(value,',',&token);CHKERRQ(ierr);
281e26ddf31SBarry Smith         ierr      = PetscTokenFind(token,&value);CHKERRQ(ierr);
282e26ddf31SBarry Smith         while (n < nmax) {
283e26ddf31SBarry Smith           if (!value) break;
284cfbddea1SSatish Balay           ierr = PetscOptionsStringToReal(value,dvalue);CHKERRQ(ierr);
285e26ddf31SBarry Smith           dvalue++;
286e26ddf31SBarry Smith           n++;
287e26ddf31SBarry Smith           ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
288e26ddf31SBarry Smith         }
2898c74ee41SBarry Smith         ierr = PetscTokenDestroy(&token);CHKERRQ(ierr);
290e26ddf31SBarry Smith       }
291e26ddf31SBarry Smith       break;
2926356e834SBarry Smith     case OPTION_INT:
293e26ddf31SBarry 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);
2943fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
2953fc1eb6aSBarry Smith       if (str[0]) {
296d25d7f95SJed Brown #if defined(PETSC_SIZEOF_LONG_LONG)
297d25d7f95SJed Brown         long long lid;
298d25d7f95SJed Brown         sscanf(str,"%lld",&lid);
299d25d7f95SJed Brown         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);
300c272547aSJed Brown #else
301d25d7f95SJed Brown         long  lid;
302d25d7f95SJed Brown         sscanf(str,"%ld",&lid);
303d25d7f95SJed Brown         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);
304c272547aSJed Brown #endif
305a297a907SKarl Rupp 
306d25d7f95SJed Brown         next->set = PETSC_TRUE;
307d25d7f95SJed Brown         *((PetscInt*)next->data) = (PetscInt)lid;
308aee2cecaSBarry Smith       }
309aee2cecaSBarry Smith       break;
310aee2cecaSBarry Smith     case OPTION_REAL:
311e26ddf31SBarry 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);
3123fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3133fc1eb6aSBarry Smith       if (str[0]) {
314ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
315a4404d99SBarry Smith         sscanf(str,"%e",&ir);
316ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL_DOUBLE)
317aee2cecaSBarry Smith         sscanf(str,"%le",&ir);
318ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL___FLOAT128)
319d9822059SBarry Smith         ir = strtoflt128(str,0);
320d9822059SBarry Smith #else
321513dbe71SLisandro Dalcin         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Unknown scalar type");
322a4404d99SBarry Smith #endif
323aee2cecaSBarry Smith         next->set                 = PETSC_TRUE;
324aee2cecaSBarry Smith         *((PetscReal*)next->data) = ir;
325aee2cecaSBarry Smith       }
326aee2cecaSBarry Smith       break;
3277781c08eSBarry Smith     case OPTION_BOOL:
3287781c08eSBarry 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);
3297781c08eSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3307781c08eSBarry Smith       if (str[0]) {
3317781c08eSBarry Smith         ierr = PetscOptionsStringToBool(str,&bid);CHKERRQ(ierr);
3327781c08eSBarry Smith         next->set = PETSC_TRUE;
3337781c08eSBarry Smith         *((PetscBool*)next->data) = bid;
3347781c08eSBarry Smith       }
3357781c08eSBarry Smith       break;
336aee2cecaSBarry Smith     case OPTION_STRING:
337e26ddf31SBarry 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);
3383fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3393fc1eb6aSBarry Smith       if (str[0]) {
340aee2cecaSBarry Smith         next->set = PETSC_TRUE;
34164facd6cSBarry Smith         /* must use system malloc since SAWs may free this */
3425b02f95dSBarry Smith         ierr = PetscStrdup(str,(char**)&next->data);CHKERRQ(ierr);
3436356e834SBarry Smith       }
3446356e834SBarry Smith       break;
345a264d7a6SBarry Smith     case OPTION_FLIST:
346140e18c1SBarry Smith       ierr = PetscFunctionListPrintTypes(PETSC_COMM_WORLD,stdout,PetscOptionsObject.prefix,next->option,next->text,next->man,next->flist,(char*)next->data);CHKERRQ(ierr);
3473cc1e11dSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3483cc1e11dSBarry Smith       if (str[0]) {
3493cc1e11dSBarry Smith         PetscOptionsObject.changedmethod = PETSC_TRUE;
3503cc1e11dSBarry Smith         next->set = PETSC_TRUE;
35164facd6cSBarry Smith         /* must use system malloc since SAWs may free this */
3525b02f95dSBarry Smith         ierr = PetscStrdup(str,(char**)&next->data);CHKERRQ(ierr);
3533cc1e11dSBarry Smith       }
3543cc1e11dSBarry Smith       break;
355b432afdaSMatthew Knepley     default:
356b432afdaSMatthew Knepley       break;
3576356e834SBarry Smith     }
3586356e834SBarry Smith     next = next->next;
3596356e834SBarry Smith   }
3606356e834SBarry Smith   PetscFunctionReturn(0);
3616356e834SBarry Smith }
3626356e834SBarry Smith 
363e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
364e04113cfSBarry Smith #include <petscviewersaws.h>
365d5649816SBarry Smith 
366d5649816SBarry Smith static int count = 0;
367d5649816SBarry Smith 
368b3506946SBarry Smith #undef __FUNCT__
369e04113cfSBarry Smith #define __FUNCT__ "PetscOptionsSAWsDestroy"
370e04113cfSBarry Smith PetscErrorCode PetscOptionsSAWsDestroy(void)
371d5649816SBarry Smith {
3722657e9d9SBarry Smith   PetscFunctionBegin;
373d5649816SBarry Smith   PetscFunctionReturn(0);
374d5649816SBarry Smith }
375d5649816SBarry Smith 
376d5649816SBarry Smith #undef __FUNCT__
3777781c08eSBarry Smith #define __FUNCT__ "PetscOptionsSAWsInput"
378b3506946SBarry Smith /*
3797781c08eSBarry Smith     PetscOptionsSAWsInput - Presents all the PETSc Options processed by the program so the user may change them at runtime using the SAWs
380b3506946SBarry Smith 
381b3506946SBarry Smith     Bugs:
382b3506946SBarry Smith +    All processes must traverse through the exact same set of option queries do to the call to PetscScanString()
383b3506946SBarry Smith .    Internal strings have arbitrary length and string copies are not checked that they fit into string space
384b3506946SBarry Smith -    Only works for PetscInt == int, PetscReal == double etc
385b3506946SBarry Smith 
386b3506946SBarry Smith 
387b3506946SBarry Smith */
388475446a1SBarry Smith PetscErrorCode PetscOptionsSAWsInput()
389b3506946SBarry Smith {
390b3506946SBarry Smith   PetscErrorCode ierr;
391b3506946SBarry Smith   PetscOptions   next     = PetscOptionsObject.next;
392d5649816SBarry Smith   static int     mancount = 0;
393b3506946SBarry Smith   char           options[16];
3947aab2a10SBarry Smith   PetscBool      changedmethod = PETSC_FALSE;
39588a9752cSBarry Smith   char           manname[16],textname[16];
3962657e9d9SBarry Smith   char           dir[1024];
397b3506946SBarry Smith 
398b3506946SBarry Smith   /* the next line is a bug, this will only work if all processors are here, the comm passed in is ignored!!! */
399b3506946SBarry Smith   sprintf(options,"Options_%d",count++);
400a297a907SKarl Rupp 
401e04113cfSBarry Smith   PetscOptionsObject.pprefix = PetscOptionsObject.prefix; /* SAWs will change this, so cannot pass prefix directly */
4021bc75a8dSBarry Smith 
403d50831c4SBarry Smith   ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s","_title");CHKERRQ(ierr);
4047781c08eSBarry Smith   PetscStackCallSAWs(SAWs_Register,(dir,&PetscOptionsObject.title,1,SAWs_READ,SAWs_STRING));
4057781c08eSBarry Smith   ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s","prefix");CHKERRQ(ierr);
4062657e9d9SBarry Smith   PetscStackCallSAWs(SAWs_Register,(dir,&PetscOptionsObject.pprefix,1,SAWs_READ,SAWs_STRING));
4072657e9d9SBarry Smith   PetscStackCallSAWs(SAWs_Register,("/PETSc/Options/ChangedMethod",&changedmethod,1,SAWs_WRITE,SAWs_BOOLEAN));
408b3506946SBarry Smith 
409b3506946SBarry Smith   while (next) {
410d50831c4SBarry Smith     sprintf(manname,"_man_%d",mancount);
4112657e9d9SBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",manname);CHKERRQ(ierr);
4127781c08eSBarry Smith     PetscStackCallSAWs(SAWs_Register,(dir,&next->man,1,SAWs_READ,SAWs_STRING));
413d50831c4SBarry Smith     sprintf(textname,"_text_%d",mancount++);
4147781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",textname);CHKERRQ(ierr);
4157781c08eSBarry Smith     PetscStackCallSAWs(SAWs_Register,(dir,&next->text,1,SAWs_READ,SAWs_STRING));
4169f32e415SBarry Smith 
417b3506946SBarry Smith     switch (next->type) {
418b3506946SBarry Smith     case OPTION_HEAD:
419b3506946SBarry Smith       break;
420b3506946SBarry Smith     case OPTION_INT_ARRAY:
4217781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4222657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_INT));
423b3506946SBarry Smith       break;
424b3506946SBarry Smith     case OPTION_REAL_ARRAY:
4257781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4262657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_DOUBLE));
427b3506946SBarry Smith       break;
428b3506946SBarry Smith     case OPTION_INT:
4297781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4302657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_INT));
431b3506946SBarry Smith       break;
432b3506946SBarry Smith     case OPTION_REAL:
4337781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4342657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_DOUBLE));
435b3506946SBarry Smith       break;
4367781c08eSBarry Smith     case OPTION_BOOL:
4377781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4382657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_BOOLEAN));
4391ae3d29cSBarry Smith       break;
4407781c08eSBarry Smith     case OPTION_BOOL_ARRAY:
4417781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4422657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_BOOLEAN));
44371f08665SBarry Smith       break;
444b3506946SBarry Smith     case OPTION_STRING:
4457781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4467781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
4471ae3d29cSBarry Smith       break;
4481ae3d29cSBarry Smith     case OPTION_STRING_ARRAY:
4497781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4502657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_STRING));
451b3506946SBarry Smith       break;
452a264d7a6SBarry Smith     case OPTION_FLIST:
453a264d7a6SBarry Smith       {
454a264d7a6SBarry Smith       PetscInt ntext;
4557781c08eSBarry Smith       ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4567781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
457a264d7a6SBarry Smith       ierr = PetscFunctionListGet(next->flist,(const char***)&next->edata,&ntext);CHKERRQ(ierr);
458a264d7a6SBarry Smith       PetscStackCallSAWs(SAWs_Set_Legal_Variable_Values,(dir,ntext,next->edata));
459a264d7a6SBarry Smith       }
460a264d7a6SBarry Smith       break;
4611ae3d29cSBarry Smith     case OPTION_ELIST:
462a264d7a6SBarry Smith       {
463a264d7a6SBarry Smith       PetscInt ntext = next->nlist;
4647781c08eSBarry Smith       ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4657781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
466315e5ce4SBarry Smith       ierr = PetscMalloc((ntext+1)*sizeof(char),&next->edata);CHKERRQ(ierr);
4671ae3d29cSBarry Smith       ierr = PetscMemcpy(next->edata,next->list,ntext*sizeof(char*));CHKERRQ(ierr);
468a264d7a6SBarry Smith       PetscStackCallSAWs(SAWs_Set_Legal_Variable_Values,(dir,ntext,next->edata));
469a264d7a6SBarry Smith       }
470a264d7a6SBarry Smith       break;
471b3506946SBarry Smith     default:
472b3506946SBarry Smith       break;
473b3506946SBarry Smith     }
474b3506946SBarry Smith     next = next->next;
475b3506946SBarry Smith   }
476b3506946SBarry Smith 
477b3506946SBarry Smith   /* wait until accessor has unlocked the memory */
4787aab2a10SBarry Smith   ierr = PetscSAWsBlock();CHKERRQ(ierr);
479b3506946SBarry Smith 
48088a9752cSBarry Smith   /* determine if any values have been set in GUI */
48188a9752cSBarry Smith   next = PetscOptionsObject.next;
48288a9752cSBarry Smith   while (next) {
48388a9752cSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
48488a9752cSBarry Smith     PetscStackCallSAWs(SAWs_Selected,(dir,&next->set));
48588a9752cSBarry Smith     next = next->next;
48688a9752cSBarry Smith   }
48788a9752cSBarry Smith 
488b3506946SBarry Smith   /* reset counter to -2; this updates the screen with the new options for the selected method */
489b3506946SBarry Smith   if (changedmethod) PetscOptionsPublishCount = -2;
490b3506946SBarry Smith 
4919a492a5cSBarry Smith   PetscStackCallSAWs(SAWs_Delete,("/PETSc/Options"));
492b3506946SBarry Smith   PetscFunctionReturn(0);
493b3506946SBarry Smith }
494b3506946SBarry Smith #endif
495b3506946SBarry Smith 
4966356e834SBarry Smith #undef __FUNCT__
49753acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEnd_Private"
49853acd3b1SBarry Smith PetscErrorCode PetscOptionsEnd_Private(void)
49953acd3b1SBarry Smith {
50053acd3b1SBarry Smith   PetscErrorCode ierr;
5016356e834SBarry Smith   PetscOptions   last;
5026356e834SBarry Smith   char           option[256],value[1024],tmp[32];
503330cf3c9SBarry Smith   size_t         j;
50453acd3b1SBarry Smith 
50553acd3b1SBarry Smith   PetscFunctionBegin;
506aee2cecaSBarry Smith   if (PetscOptionsObject.next) {
507b3506946SBarry Smith     if (!PetscOptionsPublishCount) {
508a264d7a6SBarry Smith #if defined(PETSC_HAVE_SAWS)
509475446a1SBarry Smith       ierr = PetscOptionsSAWsInput();CHKERRQ(ierr);
510b3506946SBarry Smith #else
51171f08665SBarry Smith       ierr = PetscOptionsGetFromTextInput();CHKERRQ(ierr);
512b3506946SBarry Smith #endif
513aee2cecaSBarry Smith     }
514aee2cecaSBarry Smith   }
5156356e834SBarry Smith 
516c31cb41cSBarry Smith   ierr = PetscFree(PetscOptionsObject.title);CHKERRQ(ierr);
517c31cb41cSBarry Smith   ierr = PetscFree(PetscOptionsObject.prefix);CHKERRQ(ierr);
5186356e834SBarry Smith 
5196356e834SBarry Smith   /* reset counter to -2; this updates the screen with the new options for the selected method */
5206356e834SBarry Smith   if (PetscOptionsObject.changedmethod) PetscOptionsPublishCount = -2;
52161b37b28SSatish Balay   /* reset alreadyprinted flag */
52261b37b28SSatish Balay   PetscOptionsObject.alreadyprinted = PETSC_FALSE;
5233194b578SJed Brown   if (PetscOptionsObject.object) PetscOptionsObject.object->optionsprinted = PETSC_TRUE;
5240298fd71SBarry Smith   PetscOptionsObject.object = NULL;
5256356e834SBarry Smith 
5266356e834SBarry Smith   while (PetscOptionsObject.next) {
5276356e834SBarry Smith     if (PetscOptionsObject.next->set) {
5286356e834SBarry Smith       if (PetscOptionsObject.prefix) {
5296356e834SBarry Smith         ierr = PetscStrcpy(option,"-");CHKERRQ(ierr);
5306356e834SBarry Smith         ierr = PetscStrcat(option,PetscOptionsObject.prefix);CHKERRQ(ierr);
5316356e834SBarry Smith         ierr = PetscStrcat(option,PetscOptionsObject.next->option+1);CHKERRQ(ierr);
5326356e834SBarry Smith       } else {
5336356e834SBarry Smith         ierr = PetscStrcpy(option,PetscOptionsObject.next->option);CHKERRQ(ierr);
5346356e834SBarry Smith       }
5356356e834SBarry Smith 
5366356e834SBarry Smith       switch (PetscOptionsObject.next->type) {
5376356e834SBarry Smith       case OPTION_HEAD:
5386356e834SBarry Smith         break;
539e26ddf31SBarry Smith       case OPTION_INT_ARRAY:
540e26ddf31SBarry Smith         sprintf(value,"%d",(int)((PetscInt*)PetscOptionsObject.next->data)[0]);
541e26ddf31SBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
542e26ddf31SBarry Smith           sprintf(tmp,"%d",(int)((PetscInt*)PetscOptionsObject.next->data)[j]);
543e26ddf31SBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
544e26ddf31SBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
545e26ddf31SBarry Smith         }
546e26ddf31SBarry Smith         break;
5476356e834SBarry Smith       case OPTION_INT:
5487a72a596SBarry Smith         sprintf(value,"%d",(int) *(PetscInt*)PetscOptionsObject.next->data);
5496356e834SBarry Smith         break;
5506356e834SBarry Smith       case OPTION_REAL:
5517a72a596SBarry Smith         sprintf(value,"%g",(double) *(PetscReal*)PetscOptionsObject.next->data);
5526356e834SBarry Smith         break;
5536356e834SBarry Smith       case OPTION_REAL_ARRAY:
5547a72a596SBarry Smith         sprintf(value,"%g",(double)((PetscReal*)PetscOptionsObject.next->data)[0]);
5556356e834SBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
5567a72a596SBarry Smith           sprintf(tmp,"%g",(double)((PetscReal*)PetscOptionsObject.next->data)[j]);
5576356e834SBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5586356e834SBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5596356e834SBarry Smith         }
5606356e834SBarry Smith         break;
5617781c08eSBarry Smith       case OPTION_BOOL:
56271f08665SBarry Smith         sprintf(value,"%d",*(int*)PetscOptionsObject.next->data);
5636356e834SBarry Smith         break;
5647781c08eSBarry Smith       case OPTION_BOOL_ARRAY:
565ace3abfcSBarry Smith         sprintf(value,"%d",(int)((PetscBool*)PetscOptionsObject.next->data)[0]);
5661ae3d29cSBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
567ace3abfcSBarry Smith           sprintf(tmp,"%d",(int)((PetscBool*)PetscOptionsObject.next->data)[j]);
5681ae3d29cSBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5691ae3d29cSBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5701ae3d29cSBarry Smith         }
5711ae3d29cSBarry Smith         break;
572a264d7a6SBarry Smith       case OPTION_FLIST:
5731ae3d29cSBarry Smith       case OPTION_ELIST:
5747781c08eSBarry Smith         ierr = PetscStrcpy(value,(char*)PetscOptionsObject.next->data);CHKERRQ(ierr);
5756356e834SBarry Smith         break;
5761ae3d29cSBarry Smith       case OPTION_STRING:
577475446a1SBarry Smith         ierr = PetscStrcpy(value,(char*)PetscOptionsObject.next->data);CHKERRQ(ierr);
5788da2146eSBarry Smith         break;
5791ae3d29cSBarry Smith       case OPTION_STRING_ARRAY:
5801ae3d29cSBarry Smith         sprintf(value,"%s",((char**)PetscOptionsObject.next->data)[0]);
5811ae3d29cSBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
5821ae3d29cSBarry Smith           sprintf(tmp,"%s",((char**)PetscOptionsObject.next->data)[j]);
5831ae3d29cSBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5841ae3d29cSBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5851ae3d29cSBarry Smith         }
5866356e834SBarry Smith         break;
5876356e834SBarry Smith       }
5886356e834SBarry Smith       ierr = PetscOptionsSetValue(option,value);CHKERRQ(ierr);
5896356e834SBarry Smith     }
590503cfb0cSBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->text);CHKERRQ(ierr);
591503cfb0cSBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->option);CHKERRQ(ierr);
5926356e834SBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->man);CHKERRQ(ierr);
59371f08665SBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->edata);CHKERRQ(ierr);
594c979a496SBarry Smith 
595c979a496SBarry Smith     if ((PetscOptionsObject.next->type == OPTION_STRING) || (PetscOptionsObject.next->type == OPTION_FLIST) || (PetscOptionsObject.next->type == OPTION_ELIST)){
59664facd6cSBarry Smith       /* must use system free since SAWs may have allocated it */
597c979a496SBarry Smith       free(PetscOptionsObject.next->data);
598c979a496SBarry Smith     } else {
5997781c08eSBarry Smith       ierr   = PetscFree(PetscOptionsObject.next->data);CHKERRQ(ierr);
600c979a496SBarry Smith     }
6017781c08eSBarry Smith 
6026356e834SBarry Smith     last                    = PetscOptionsObject.next;
6036356e834SBarry Smith     PetscOptionsObject.next = PetscOptionsObject.next->next;
6046356e834SBarry Smith     ierr                    = PetscFree(last);CHKERRQ(ierr);
6056356e834SBarry Smith   }
6066356e834SBarry Smith   PetscOptionsObject.next = 0;
60753acd3b1SBarry Smith   PetscFunctionReturn(0);
60853acd3b1SBarry Smith }
60953acd3b1SBarry Smith 
61053acd3b1SBarry Smith #undef __FUNCT__
61153acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEnum"
61253acd3b1SBarry Smith /*@C
61353acd3b1SBarry Smith    PetscOptionsEnum - Gets the enum value for a particular option in the database.
61453acd3b1SBarry Smith 
6153f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
61653acd3b1SBarry Smith 
61753acd3b1SBarry Smith    Input Parameters:
61853acd3b1SBarry Smith +  opt - option name
61953acd3b1SBarry Smith .  text - short string that describes the option
62053acd3b1SBarry Smith .  man - manual page with additional information on option
62153acd3b1SBarry Smith .  list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null
622*0fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
623*0fdccdaeSBarry Smith $                 PetscOptionsEnum(..., obj->value,&object->value,...) or
624*0fdccdaeSBarry Smith $                 value = defaultvalue
625*0fdccdaeSBarry Smith $                 PetscOptionsEnum(..., value,&value,&flg);
626*0fdccdaeSBarry Smith $                 if (flg) {
62753acd3b1SBarry Smith 
62853acd3b1SBarry Smith    Output Parameter:
62953acd3b1SBarry Smith +  value - the  value to return
630b32e0204SMatthew G Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
63153acd3b1SBarry Smith 
63253acd3b1SBarry Smith    Level: beginner
63353acd3b1SBarry Smith 
63453acd3b1SBarry Smith    Concepts: options database
63553acd3b1SBarry Smith 
63653acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
63753acd3b1SBarry Smith 
63853acd3b1SBarry Smith           list is usually something like PCASMTypes or some other predefined list of enum names
63953acd3b1SBarry Smith 
64053acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
641acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
642acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
64353acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
64453acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
645acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
646a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
64753acd3b1SBarry Smith @*/
648*0fdccdaeSBarry Smith PetscErrorCode  PetscOptionsEnum(const char opt[],const char text[],const char man[],const char *const *list,PetscEnum currentvalue,PetscEnum *value,PetscBool  *set)
64953acd3b1SBarry Smith {
65053acd3b1SBarry Smith   PetscErrorCode ierr;
65153acd3b1SBarry Smith   PetscInt       ntext = 0;
652aa5bb8c0SSatish Balay   PetscInt       tval;
653ace3abfcSBarry Smith   PetscBool      tflg;
65453acd3b1SBarry Smith 
65553acd3b1SBarry Smith   PetscFunctionBegin;
65653acd3b1SBarry Smith   while (list[ntext++]) {
657e32f2f54SBarry Smith     if (ntext > 50) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument appears to be wrong or have more than 50 entries");
65853acd3b1SBarry Smith   }
659e32f2f54SBarry Smith   if (ntext < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument must have at least two entries: typename and type prefix");
66053acd3b1SBarry Smith   ntext -= 3;
661*0fdccdaeSBarry Smith   ierr   = PetscOptionsEList(opt,text,man,list,ntext,list[currentvalue],&tval,&tflg);CHKERRQ(ierr);
662aa5bb8c0SSatish Balay   /* with PETSC_USE_64BIT_INDICES sizeof(PetscInt) != sizeof(PetscEnum) */
663aa5bb8c0SSatish Balay   if (tflg) *value = (PetscEnum)tval;
664aa5bb8c0SSatish Balay   if (set)  *set   = tflg;
66553acd3b1SBarry Smith   PetscFunctionReturn(0);
66653acd3b1SBarry Smith }
66753acd3b1SBarry Smith 
66853acd3b1SBarry Smith /* -------------------------------------------------------------------------------------------------------------*/
66953acd3b1SBarry Smith #undef __FUNCT__
67053acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsInt"
67153acd3b1SBarry Smith /*@C
67253acd3b1SBarry Smith    PetscOptionsInt - Gets the integer value for a particular option in the database.
67353acd3b1SBarry Smith 
6743f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
67553acd3b1SBarry Smith 
67653acd3b1SBarry Smith    Input Parameters:
67753acd3b1SBarry Smith +  opt - option name
67853acd3b1SBarry Smith .  text - short string that describes the option
67953acd3b1SBarry Smith .  man - manual page with additional information on option
680*0fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
681*0fdccdaeSBarry Smith $                 PetscOptionsInt(..., obj->value,&object->value,...) or
682*0fdccdaeSBarry Smith $                 value = defaultvalue
683*0fdccdaeSBarry Smith $                 PetscOptionsInt(..., value,&value,&flg);
684*0fdccdaeSBarry Smith $                 if (flg) {
68553acd3b1SBarry Smith 
68653acd3b1SBarry Smith    Output Parameter:
68753acd3b1SBarry Smith +  value - the integer value to return
68853acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
68953acd3b1SBarry Smith 
69053acd3b1SBarry Smith    Level: beginner
69153acd3b1SBarry Smith 
69253acd3b1SBarry Smith    Concepts: options database^has int
69353acd3b1SBarry Smith 
69453acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
69553acd3b1SBarry Smith 
69653acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
697acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
698acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
69953acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
70053acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
701acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
702a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
70353acd3b1SBarry Smith @*/
704*0fdccdaeSBarry Smith PetscErrorCode  PetscOptionsInt(const char opt[],const char text[],const char man[],PetscInt currentvalue,PetscInt *value,PetscBool  *set)
70553acd3b1SBarry Smith {
70653acd3b1SBarry Smith   PetscErrorCode ierr;
7076356e834SBarry Smith   PetscOptions   amsopt;
70853acd3b1SBarry Smith 
70953acd3b1SBarry Smith   PetscFunctionBegin;
710b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
7116356e834SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_INT,&amsopt);CHKERRQ(ierr);
7126356e834SBarry Smith     ierr = PetscMalloc(sizeof(PetscInt),&amsopt->data);CHKERRQ(ierr);
713a297a907SKarl Rupp 
714*0fdccdaeSBarry Smith     *(PetscInt*)amsopt->data = currentvalue;
715af6d86caSBarry Smith   }
71653acd3b1SBarry Smith   ierr = PetscOptionsGetInt(PetscOptionsObject.prefix,opt,value,set);CHKERRQ(ierr);
71761b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
718*0fdccdaeSBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%d>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,currentvalue,text,ManSection(man));CHKERRQ(ierr);
71953acd3b1SBarry Smith   }
72053acd3b1SBarry Smith   PetscFunctionReturn(0);
72153acd3b1SBarry Smith }
72253acd3b1SBarry Smith 
72353acd3b1SBarry Smith #undef __FUNCT__
72453acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsString"
72553acd3b1SBarry Smith /*@C
72653acd3b1SBarry Smith    PetscOptionsString - Gets the string value for a particular option in the database.
72753acd3b1SBarry Smith 
7283f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
72953acd3b1SBarry Smith 
73053acd3b1SBarry Smith    Input Parameters:
73153acd3b1SBarry Smith +  opt - option name
73253acd3b1SBarry Smith .  text - short string that describes the option
73353acd3b1SBarry Smith .  man - manual page with additional information on option
734*0fdccdaeSBarry Smith .  currentvalue - the current value; caller is responsible for setting this value correctly. This is not used to set value
735bcbf2dc5SJed Brown -  len - length of the result string including null terminator
73653acd3b1SBarry Smith 
73753acd3b1SBarry Smith    Output Parameter:
73853acd3b1SBarry Smith +  value - the value to return
73953acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
74053acd3b1SBarry Smith 
74153acd3b1SBarry Smith    Level: beginner
74253acd3b1SBarry Smith 
74353acd3b1SBarry Smith    Concepts: options database^has int
74453acd3b1SBarry Smith 
74553acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
74653acd3b1SBarry Smith 
7477fccdfe4SBarry 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).
7487fccdfe4SBarry Smith 
74953acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
750acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
751acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsReal(), PetscOptionsBool(),
75253acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
75353acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
754acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
755a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
75653acd3b1SBarry Smith @*/
757*0fdccdaeSBarry Smith PetscErrorCode  PetscOptionsString(const char opt[],const char text[],const char man[],const char currentvalue[],char value[],size_t len,PetscBool  *set)
75853acd3b1SBarry Smith {
75953acd3b1SBarry Smith   PetscErrorCode ierr;
760aee2cecaSBarry Smith   PetscOptions   amsopt;
76153acd3b1SBarry Smith 
76253acd3b1SBarry Smith   PetscFunctionBegin;
763b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
764aee2cecaSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
76564facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
766*0fdccdaeSBarry Smith     ierr = PetscStrdup(currentvalue ? currentvalue : "",(char**)&amsopt->data);CHKERRQ(ierr);
767af6d86caSBarry Smith   }
76853acd3b1SBarry Smith   ierr = PetscOptionsGetString(PetscOptionsObject.prefix,opt,value,len,set);CHKERRQ(ierr);
76961b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
770*0fdccdaeSBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%s>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,currentvalue,text,ManSection(man));CHKERRQ(ierr);
77153acd3b1SBarry Smith   }
77253acd3b1SBarry Smith   PetscFunctionReturn(0);
77353acd3b1SBarry Smith }
77453acd3b1SBarry Smith 
77553acd3b1SBarry Smith #undef __FUNCT__
77653acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsReal"
77753acd3b1SBarry Smith /*@C
77853acd3b1SBarry Smith    PetscOptionsReal - Gets the PetscReal value for a particular option in the database.
77953acd3b1SBarry Smith 
7803f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
78153acd3b1SBarry Smith 
78253acd3b1SBarry Smith    Input Parameters:
78353acd3b1SBarry Smith +  opt - option name
78453acd3b1SBarry Smith .  text - short string that describes the option
78553acd3b1SBarry Smith .  man - manual page with additional information on option
786*0fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
787*0fdccdaeSBarry Smith $                 PetscOptionsReal(..., obj->value,&object->value,...) or
788*0fdccdaeSBarry Smith $                 value = defaultvalue
789*0fdccdaeSBarry Smith $                 PetscOptionsReal(..., value,&value,&flg);
790*0fdccdaeSBarry Smith $                 if (flg) {
79153acd3b1SBarry Smith 
79253acd3b1SBarry Smith    Output Parameter:
79353acd3b1SBarry Smith +  value - the value to return
79453acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
79553acd3b1SBarry Smith 
79653acd3b1SBarry Smith    Level: beginner
79753acd3b1SBarry Smith 
79853acd3b1SBarry Smith    Concepts: options database^has int
79953acd3b1SBarry Smith 
80053acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
80153acd3b1SBarry Smith 
80253acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
803acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
804acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
80553acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
80653acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
807acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
808a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
80953acd3b1SBarry Smith @*/
810*0fdccdaeSBarry Smith PetscErrorCode  PetscOptionsReal(const char opt[],const char text[],const char man[],PetscReal currentvalue,PetscReal *value,PetscBool  *set)
81153acd3b1SBarry Smith {
81253acd3b1SBarry Smith   PetscErrorCode ierr;
813538aa990SBarry Smith   PetscOptions   amsopt;
81453acd3b1SBarry Smith 
81553acd3b1SBarry Smith   PetscFunctionBegin;
816b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
817538aa990SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_REAL,&amsopt);CHKERRQ(ierr);
818538aa990SBarry Smith     ierr = PetscMalloc(sizeof(PetscReal),&amsopt->data);CHKERRQ(ierr);
819a297a907SKarl Rupp 
820*0fdccdaeSBarry Smith     *(PetscReal*)amsopt->data = currentvalue;
821538aa990SBarry Smith   }
82253acd3b1SBarry Smith   ierr = PetscOptionsGetReal(PetscOptionsObject.prefix,opt,value,set);CHKERRQ(ierr);
82361b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
824*0fdccdaeSBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%g>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,(double)currentvalue,text,ManSection(man));CHKERRQ(ierr);
82553acd3b1SBarry Smith   }
82653acd3b1SBarry Smith   PetscFunctionReturn(0);
82753acd3b1SBarry Smith }
82853acd3b1SBarry Smith 
82953acd3b1SBarry Smith #undef __FUNCT__
83053acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsScalar"
83153acd3b1SBarry Smith /*@C
83253acd3b1SBarry Smith    PetscOptionsScalar - Gets the scalar value for a particular option in the database.
83353acd3b1SBarry Smith 
8343f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
83553acd3b1SBarry Smith 
83653acd3b1SBarry Smith    Input Parameters:
83753acd3b1SBarry Smith +  opt - option name
83853acd3b1SBarry Smith .  text - short string that describes the option
83953acd3b1SBarry Smith .  man - manual page with additional information on option
840*0fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
841*0fdccdaeSBarry Smith $                 PetscOptionsScalar(..., obj->value,&object->value,...) or
842*0fdccdaeSBarry Smith $                 value = defaultvalue
843*0fdccdaeSBarry Smith $                 PetscOptionsScalar(..., value,&value,&flg);
844*0fdccdaeSBarry Smith $                 if (flg) {
845*0fdccdaeSBarry Smith 
84653acd3b1SBarry Smith 
84753acd3b1SBarry Smith    Output Parameter:
84853acd3b1SBarry Smith +  value - the value to return
84953acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
85053acd3b1SBarry Smith 
85153acd3b1SBarry Smith    Level: beginner
85253acd3b1SBarry Smith 
85353acd3b1SBarry Smith    Concepts: options database^has int
85453acd3b1SBarry Smith 
85553acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
85653acd3b1SBarry Smith 
85753acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
858acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
859acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
86053acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
86153acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
862acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
863a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
86453acd3b1SBarry Smith @*/
865*0fdccdaeSBarry Smith PetscErrorCode  PetscOptionsScalar(const char opt[],const char text[],const char man[],PetscScalar currentvalue,PetscScalar *value,PetscBool  *set)
86653acd3b1SBarry Smith {
86753acd3b1SBarry Smith   PetscErrorCode ierr;
86853acd3b1SBarry Smith 
86953acd3b1SBarry Smith   PetscFunctionBegin;
87053acd3b1SBarry Smith #if !defined(PETSC_USE_COMPLEX)
871*0fdccdaeSBarry Smith   ierr = PetscOptionsReal(opt,text,man,currentvalue,value,set);CHKERRQ(ierr);
87253acd3b1SBarry Smith #else
87353acd3b1SBarry Smith   ierr = PetscOptionsGetScalar(PetscOptionsObject.prefix,opt,value,set);CHKERRQ(ierr);
87453acd3b1SBarry Smith #endif
87553acd3b1SBarry Smith   PetscFunctionReturn(0);
87653acd3b1SBarry Smith }
87753acd3b1SBarry Smith 
87853acd3b1SBarry Smith #undef __FUNCT__
87953acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsName"
88053acd3b1SBarry Smith /*@C
88190d69ab7SBarry 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
88290d69ab7SBarry Smith                       its value is set to false.
88353acd3b1SBarry Smith 
8843f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
88553acd3b1SBarry Smith 
88653acd3b1SBarry Smith    Input Parameters:
88753acd3b1SBarry Smith +  opt - option name
88853acd3b1SBarry Smith .  text - short string that describes the option
88953acd3b1SBarry Smith -  man - manual page with additional information on option
89053acd3b1SBarry Smith 
89153acd3b1SBarry Smith    Output Parameter:
89253acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
89353acd3b1SBarry Smith 
89453acd3b1SBarry Smith    Level: beginner
89553acd3b1SBarry Smith 
89653acd3b1SBarry Smith    Concepts: options database^has int
89753acd3b1SBarry Smith 
89853acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
89953acd3b1SBarry Smith 
90053acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
901acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
902acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
90353acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
90453acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
905acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
906a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
90753acd3b1SBarry Smith @*/
9087087cfbeSBarry Smith PetscErrorCode  PetscOptionsName(const char opt[],const char text[],const char man[],PetscBool  *flg)
90953acd3b1SBarry Smith {
91053acd3b1SBarry Smith   PetscErrorCode ierr;
9111ae3d29cSBarry Smith   PetscOptions   amsopt;
91253acd3b1SBarry Smith 
91353acd3b1SBarry Smith   PetscFunctionBegin;
9141ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
9157781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
916ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
917a297a907SKarl Rupp 
918ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
9191ae3d29cSBarry Smith   }
92053acd3b1SBarry Smith   ierr = PetscOptionsHasName(PetscOptionsObject.prefix,opt,flg);CHKERRQ(ierr);
92161b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
9222aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
92353acd3b1SBarry Smith   }
92453acd3b1SBarry Smith   PetscFunctionReturn(0);
92553acd3b1SBarry Smith }
92653acd3b1SBarry Smith 
92753acd3b1SBarry Smith #undef __FUNCT__
928a264d7a6SBarry Smith #define __FUNCT__ "PetscOptionsFList"
92953acd3b1SBarry Smith /*@C
930a264d7a6SBarry Smith      PetscOptionsFList - Puts a list of option values that a single one may be selected from
93153acd3b1SBarry Smith 
9323f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
93353acd3b1SBarry Smith 
93453acd3b1SBarry Smith    Input Parameters:
93553acd3b1SBarry Smith +  opt - option name
93653acd3b1SBarry Smith .  text - short string that describes the option
93753acd3b1SBarry Smith .  man - manual page with additional information on option
93853acd3b1SBarry Smith .  list - the possible choices
939*0fdccdaeSBarry Smith .  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with
940*0fdccdaeSBarry Smith $                 PetscOptionsFlist(..., obj->value,value,len,&flg);
941*0fdccdaeSBarry Smith $                 if (flg) {
9423cc1e11dSBarry Smith -  len - the length of the character array value
94353acd3b1SBarry Smith 
94453acd3b1SBarry Smith    Output Parameter:
94553acd3b1SBarry Smith +  value - the value to return
94653acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
94753acd3b1SBarry Smith 
94853acd3b1SBarry Smith    Level: intermediate
94953acd3b1SBarry Smith 
95053acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
95153acd3b1SBarry Smith 
95253acd3b1SBarry Smith    See PetscOptionsEList() for when the choices are given in a string array
95353acd3b1SBarry Smith 
95453acd3b1SBarry Smith    To get a listing of all currently specified options,
95588c29154SBarry Smith     see PetscOptionsView() or PetscOptionsGetAll()
95653acd3b1SBarry Smith 
957eabe10d7SBarry Smith    Developer Note: This cannot check for invalid selection because of things like MATAIJ that are not included in the list
958eabe10d7SBarry Smith 
95953acd3b1SBarry Smith    Concepts: options database^list
96053acd3b1SBarry Smith 
96153acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
962acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
96353acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
96453acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
965acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
966a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsEnum()
96753acd3b1SBarry Smith @*/
968*0fdccdaeSBarry Smith PetscErrorCode  PetscOptionsFList(const char opt[],const char ltext[],const char man[],PetscFunctionList list,const char currentvalue[],char value[],size_t len,PetscBool  *set)
96953acd3b1SBarry Smith {
97053acd3b1SBarry Smith   PetscErrorCode ierr;
9713cc1e11dSBarry Smith   PetscOptions   amsopt;
97253acd3b1SBarry Smith 
97353acd3b1SBarry Smith   PetscFunctionBegin;
974b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
975a264d7a6SBarry Smith     ierr = PetscOptionsCreate_Private(opt,ltext,man,OPTION_FLIST,&amsopt);CHKERRQ(ierr);
97664facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
977*0fdccdaeSBarry Smith     ierr = PetscStrdup(currentvalue ? currentvalue : "",(char**)&amsopt->data);CHKERRQ(ierr);
9783cc1e11dSBarry Smith     amsopt->flist = list;
9793cc1e11dSBarry Smith   }
98053acd3b1SBarry Smith   ierr = PetscOptionsGetString(PetscOptionsObject.prefix,opt,value,len,set);CHKERRQ(ierr);
98161b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
982*0fdccdaeSBarry Smith     ierr = PetscFunctionListPrintTypes(PetscOptionsObject.comm,stdout,PetscOptionsObject.prefix,opt,ltext,man,list,currentvalue);CHKERRQ(ierr);CHKERRQ(ierr);
98353acd3b1SBarry Smith   }
98453acd3b1SBarry Smith   PetscFunctionReturn(0);
98553acd3b1SBarry Smith }
98653acd3b1SBarry Smith 
98753acd3b1SBarry Smith #undef __FUNCT__
98853acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEList"
98953acd3b1SBarry Smith /*@C
99053acd3b1SBarry Smith      PetscOptionsEList - Puts a list of option values that a single one may be selected from
99153acd3b1SBarry Smith 
9923f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
99353acd3b1SBarry Smith 
99453acd3b1SBarry Smith    Input Parameters:
99553acd3b1SBarry Smith +  opt - option name
99653acd3b1SBarry Smith .  ltext - short string that describes the option
99753acd3b1SBarry Smith .  man - manual page with additional information on option
998a264d7a6SBarry Smith .  list - the possible choices (one of these must be selected, anything else is invalid)
99953acd3b1SBarry Smith .  ntext - number of choices
1000*0fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with
1001*0fdccdaeSBarry Smith $                 PetscOptionsElist(..., obj->value,&value,&flg);
1002*0fdccdaeSBarry Smith $                 if (flg) {
1003*0fdccdaeSBarry Smith 
100453acd3b1SBarry Smith 
100553acd3b1SBarry Smith    Output Parameter:
100653acd3b1SBarry Smith +  value - the index of the value to return
100753acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
100853acd3b1SBarry Smith 
100953acd3b1SBarry Smith    Level: intermediate
101053acd3b1SBarry Smith 
101153acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
101253acd3b1SBarry Smith 
1013a264d7a6SBarry Smith    See PetscOptionsFList() for when the choices are given in a PetscFunctionList()
101453acd3b1SBarry Smith 
101553acd3b1SBarry Smith    Concepts: options database^list
101653acd3b1SBarry Smith 
101753acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1018acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
101953acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
102053acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1021acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1022a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEnum()
102353acd3b1SBarry Smith @*/
1024*0fdccdaeSBarry Smith PetscErrorCode  PetscOptionsEList(const char opt[],const char ltext[],const char man[],const char *const *list,PetscInt ntext,const char currentvalue[],PetscInt *value,PetscBool  *set)
102553acd3b1SBarry Smith {
102653acd3b1SBarry Smith   PetscErrorCode ierr;
102753acd3b1SBarry Smith   PetscInt       i;
10281ae3d29cSBarry Smith   PetscOptions   amsopt;
102953acd3b1SBarry Smith 
103053acd3b1SBarry Smith   PetscFunctionBegin;
10311ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
1032d5649816SBarry Smith     ierr = PetscOptionsCreate_Private(opt,ltext,man,OPTION_ELIST,&amsopt);CHKERRQ(ierr);
103364facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
1034*0fdccdaeSBarry Smith     ierr = PetscStrdup(currentvalue ? currentvalue : "",(char**)&amsopt->data);CHKERRQ(ierr);
10351ae3d29cSBarry Smith     amsopt->list  = list;
10361ae3d29cSBarry Smith     amsopt->nlist = ntext;
10371ae3d29cSBarry Smith   }
103853acd3b1SBarry Smith   ierr = PetscOptionsGetEList(PetscOptionsObject.prefix,opt,list,ntext,value,set);CHKERRQ(ierr);
103961b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
1040*0fdccdaeSBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%s> (choose one of)",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,currentvalue);CHKERRQ(ierr);
104153acd3b1SBarry Smith     for (i=0; i<ntext; i++) {
104253acd3b1SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm," %s",list[i]);CHKERRQ(ierr);
104353acd3b1SBarry Smith     }
10442aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm," (%s)\n",ManSection(man));CHKERRQ(ierr);
104553acd3b1SBarry Smith   }
104653acd3b1SBarry Smith   PetscFunctionReturn(0);
104753acd3b1SBarry Smith }
104853acd3b1SBarry Smith 
104953acd3b1SBarry Smith #undef __FUNCT__
1050acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolGroupBegin"
105153acd3b1SBarry Smith /*@C
1052acfcf0e5SJed Brown      PetscOptionsBoolGroupBegin - First in a series of logical queries on the options database for
1053d5649816SBarry Smith        which at most a single value can be true.
105453acd3b1SBarry Smith 
10553f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
105653acd3b1SBarry Smith 
105753acd3b1SBarry Smith    Input Parameters:
105853acd3b1SBarry Smith +  opt - option name
105953acd3b1SBarry Smith .  text - short string that describes the option
106053acd3b1SBarry Smith -  man - manual page with additional information on option
106153acd3b1SBarry Smith 
106253acd3b1SBarry Smith    Output Parameter:
106353acd3b1SBarry Smith .  flg - whether that option was set or not
106453acd3b1SBarry Smith 
106553acd3b1SBarry Smith    Level: intermediate
106653acd3b1SBarry Smith 
106753acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
106853acd3b1SBarry Smith 
1069acfcf0e5SJed Brown    Must be followed by 0 or more PetscOptionsBoolGroup()s and PetscOptionsBoolGroupEnd()
107053acd3b1SBarry Smith 
107153acd3b1SBarry Smith     Concepts: options database^logical group
107253acd3b1SBarry Smith 
107353acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1074acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
107553acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
107653acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1077acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1078a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
107953acd3b1SBarry Smith @*/
10807087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolGroupBegin(const char opt[],const char text[],const char man[],PetscBool  *flg)
108153acd3b1SBarry Smith {
108253acd3b1SBarry Smith   PetscErrorCode ierr;
10831ae3d29cSBarry Smith   PetscOptions   amsopt;
108453acd3b1SBarry Smith 
108553acd3b1SBarry Smith   PetscFunctionBegin;
10861ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
10877781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1088ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1089a297a907SKarl Rupp 
1090ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
10911ae3d29cSBarry Smith   }
109268b16fdaSBarry Smith   *flg = PETSC_FALSE;
10930298fd71SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,NULL);CHKERRQ(ierr);
109461b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
109553acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  Pick at most one of -------------\n");CHKERRQ(ierr);
10962aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"    -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
109753acd3b1SBarry Smith   }
109853acd3b1SBarry Smith   PetscFunctionReturn(0);
109953acd3b1SBarry Smith }
110053acd3b1SBarry Smith 
110153acd3b1SBarry Smith #undef __FUNCT__
1102acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolGroup"
110353acd3b1SBarry Smith /*@C
1104acfcf0e5SJed Brown      PetscOptionsBoolGroup - One in a series of logical queries on the options database for
1105d5649816SBarry Smith        which at most a single value can be true.
110653acd3b1SBarry Smith 
11073f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
110853acd3b1SBarry Smith 
110953acd3b1SBarry Smith    Input Parameters:
111053acd3b1SBarry Smith +  opt - option name
111153acd3b1SBarry Smith .  text - short string that describes the option
111253acd3b1SBarry Smith -  man - manual page with additional information on option
111353acd3b1SBarry Smith 
111453acd3b1SBarry Smith    Output Parameter:
111553acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
111653acd3b1SBarry Smith 
111753acd3b1SBarry Smith    Level: intermediate
111853acd3b1SBarry Smith 
111953acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
112053acd3b1SBarry Smith 
1121acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin() and preceded a PetscOptionsBoolGroupEnd()
112253acd3b1SBarry Smith 
112353acd3b1SBarry Smith     Concepts: options database^logical group
112453acd3b1SBarry Smith 
112553acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1126acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
112753acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
112853acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1129acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1130a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
113153acd3b1SBarry Smith @*/
11327087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolGroup(const char opt[],const char text[],const char man[],PetscBool  *flg)
113353acd3b1SBarry Smith {
113453acd3b1SBarry Smith   PetscErrorCode ierr;
11351ae3d29cSBarry Smith   PetscOptions   amsopt;
113653acd3b1SBarry Smith 
113753acd3b1SBarry Smith   PetscFunctionBegin;
11381ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
11397781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1140ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1141a297a907SKarl Rupp 
1142ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
11431ae3d29cSBarry Smith   }
114417326d04SJed Brown   *flg = PETSC_FALSE;
11450298fd71SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,NULL);CHKERRQ(ierr);
114661b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
11472aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"    -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
114853acd3b1SBarry Smith   }
114953acd3b1SBarry Smith   PetscFunctionReturn(0);
115053acd3b1SBarry Smith }
115153acd3b1SBarry Smith 
115253acd3b1SBarry Smith #undef __FUNCT__
1153acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolGroupEnd"
115453acd3b1SBarry Smith /*@C
1155acfcf0e5SJed Brown      PetscOptionsBoolGroupEnd - Last in a series of logical queries on the options database for
1156d5649816SBarry Smith        which at most a single value can be true.
115753acd3b1SBarry Smith 
11583f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
115953acd3b1SBarry Smith 
116053acd3b1SBarry Smith    Input Parameters:
116153acd3b1SBarry Smith +  opt - option name
116253acd3b1SBarry Smith .  text - short string that describes the option
116353acd3b1SBarry Smith -  man - manual page with additional information on option
116453acd3b1SBarry Smith 
116553acd3b1SBarry Smith    Output Parameter:
116653acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
116753acd3b1SBarry Smith 
116853acd3b1SBarry Smith    Level: intermediate
116953acd3b1SBarry Smith 
117053acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
117153acd3b1SBarry Smith 
1172acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin()
117353acd3b1SBarry Smith 
117453acd3b1SBarry Smith     Concepts: options database^logical group
117553acd3b1SBarry Smith 
117653acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1177acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
117853acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
117953acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1180acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1181a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
118253acd3b1SBarry Smith @*/
11837087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolGroupEnd(const char opt[],const char text[],const char man[],PetscBool  *flg)
118453acd3b1SBarry Smith {
118553acd3b1SBarry Smith   PetscErrorCode ierr;
11861ae3d29cSBarry Smith   PetscOptions   amsopt;
118753acd3b1SBarry Smith 
118853acd3b1SBarry Smith   PetscFunctionBegin;
11891ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
11907781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1191ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1192a297a907SKarl Rupp 
1193ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
11941ae3d29cSBarry Smith   }
119517326d04SJed Brown   *flg = PETSC_FALSE;
11960298fd71SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,NULL);CHKERRQ(ierr);
119761b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
11982aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"    -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
119953acd3b1SBarry Smith   }
120053acd3b1SBarry Smith   PetscFunctionReturn(0);
120153acd3b1SBarry Smith }
120253acd3b1SBarry Smith 
120353acd3b1SBarry Smith #undef __FUNCT__
1204acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBool"
120553acd3b1SBarry Smith /*@C
1206acfcf0e5SJed Brown    PetscOptionsBool - Determines if a particular option is in the database with a true or false
120753acd3b1SBarry Smith 
12083f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
120953acd3b1SBarry Smith 
121053acd3b1SBarry Smith    Input Parameters:
121153acd3b1SBarry Smith +  opt - option name
121253acd3b1SBarry Smith .  text - short string that describes the option
1213868c398cSBarry Smith .  man - manual page with additional information on option
1214868c398cSBarry Smith -  deflt - the default value, if the user does not set a value then this value is returned in flg
121553acd3b1SBarry Smith 
121653acd3b1SBarry Smith    Output Parameter:
121753acd3b1SBarry Smith .  flg - PETSC_TRUE or PETSC_FALSE
121853acd3b1SBarry Smith .  set - PETSC_TRUE if found, else PETSC_FALSE
121953acd3b1SBarry Smith 
122053acd3b1SBarry Smith    Level: beginner
122153acd3b1SBarry Smith 
122253acd3b1SBarry Smith    Concepts: options database^logical
122353acd3b1SBarry Smith 
122453acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
122553acd3b1SBarry Smith 
122653acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
1227acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
1228acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
122953acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
123053acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1231acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1232a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
123353acd3b1SBarry Smith @*/
12347087cfbeSBarry Smith PetscErrorCode  PetscOptionsBool(const char opt[],const char text[],const char man[],PetscBool deflt,PetscBool  *flg,PetscBool  *set)
123553acd3b1SBarry Smith {
123653acd3b1SBarry Smith   PetscErrorCode ierr;
1237ace3abfcSBarry Smith   PetscBool      iset;
1238aee2cecaSBarry Smith   PetscOptions   amsopt;
123953acd3b1SBarry Smith 
124053acd3b1SBarry Smith   PetscFunctionBegin;
1241b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
12427781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1243ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1244a297a907SKarl Rupp 
1245ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = deflt;
1246af6d86caSBarry Smith   }
1247acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,&iset);CHKERRQ(ierr);
124853acd3b1SBarry Smith   if (!iset) {
124953acd3b1SBarry Smith     if (flg) *flg = deflt;
125053acd3b1SBarry Smith   }
125153acd3b1SBarry Smith   if (set) *set = iset;
125261b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
1253ace3abfcSBarry Smith     const char *v = PetscBools[deflt];
12542aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s: <%s> %s (%s)\n",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,v,text,ManSection(man));CHKERRQ(ierr);
125553acd3b1SBarry Smith   }
125653acd3b1SBarry Smith   PetscFunctionReturn(0);
125753acd3b1SBarry Smith }
125853acd3b1SBarry Smith 
125953acd3b1SBarry Smith #undef __FUNCT__
126053acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsRealArray"
126153acd3b1SBarry Smith /*@C
126253acd3b1SBarry Smith    PetscOptionsRealArray - Gets an array of double values for a particular
126353acd3b1SBarry Smith    option in the database. The values must be separated with commas with
126453acd3b1SBarry Smith    no intervening spaces.
126553acd3b1SBarry Smith 
12663f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
126753acd3b1SBarry Smith 
126853acd3b1SBarry Smith    Input Parameters:
126953acd3b1SBarry Smith +  opt - the option one is seeking
127053acd3b1SBarry Smith .  text - short string describing option
127153acd3b1SBarry Smith .  man - manual page for option
127253acd3b1SBarry Smith -  nmax - maximum number of values
127353acd3b1SBarry Smith 
127453acd3b1SBarry Smith    Output Parameter:
127553acd3b1SBarry Smith +  value - location to copy values
127653acd3b1SBarry Smith .  nmax - actual number of values found
127753acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
127853acd3b1SBarry Smith 
127953acd3b1SBarry Smith    Level: beginner
128053acd3b1SBarry Smith 
128153acd3b1SBarry Smith    Notes:
128253acd3b1SBarry Smith    The user should pass in an array of doubles
128353acd3b1SBarry Smith 
128453acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
128553acd3b1SBarry Smith 
128653acd3b1SBarry Smith    Concepts: options database^array of strings
128753acd3b1SBarry Smith 
128853acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1289acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
129053acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
129153acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1292acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1293a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
129453acd3b1SBarry Smith @*/
12957087cfbeSBarry Smith PetscErrorCode  PetscOptionsRealArray(const char opt[],const char text[],const char man[],PetscReal value[],PetscInt *n,PetscBool  *set)
129653acd3b1SBarry Smith {
129753acd3b1SBarry Smith   PetscErrorCode ierr;
129853acd3b1SBarry Smith   PetscInt       i;
1299e26ddf31SBarry Smith   PetscOptions   amsopt;
130053acd3b1SBarry Smith 
130153acd3b1SBarry Smith   PetscFunctionBegin;
1302b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
1303e26ddf31SBarry Smith     PetscReal *vals;
1304e26ddf31SBarry Smith 
1305e26ddf31SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_REAL_ARRAY,&amsopt);CHKERRQ(ierr);
13066e02237eSPeter Brune     ierr = PetscMalloc1((*n),(PetscReal**)&amsopt->data);CHKERRQ(ierr);
1307e26ddf31SBarry Smith     vals = (PetscReal*)amsopt->data;
1308e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1309e26ddf31SBarry Smith     amsopt->arraylength = *n;
1310e26ddf31SBarry Smith   }
131153acd3b1SBarry Smith   ierr = PetscOptionsGetRealArray(PetscOptionsObject.prefix,opt,value,n,set);CHKERRQ(ierr);
131261b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
131357622a8eSBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%g",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,(double)value[0]);CHKERRQ(ierr);
131453acd3b1SBarry Smith     for (i=1; i<*n; i++) {
131557622a8eSBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,",%g",(double)value[i]);CHKERRQ(ierr);
131653acd3b1SBarry Smith     }
13172aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
131853acd3b1SBarry Smith   }
131953acd3b1SBarry Smith   PetscFunctionReturn(0);
132053acd3b1SBarry Smith }
132153acd3b1SBarry Smith 
132253acd3b1SBarry Smith 
132353acd3b1SBarry Smith #undef __FUNCT__
132453acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsIntArray"
132553acd3b1SBarry Smith /*@C
132653acd3b1SBarry Smith    PetscOptionsIntArray - Gets an array of integers for a particular
1327b32a342fSShri Abhyankar    option in the database.
132853acd3b1SBarry Smith 
13293f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
133053acd3b1SBarry Smith 
133153acd3b1SBarry Smith    Input Parameters:
133253acd3b1SBarry Smith +  opt - the option one is seeking
133353acd3b1SBarry Smith .  text - short string describing option
133453acd3b1SBarry Smith .  man - manual page for option
1335f8a50e2bSBarry Smith -  n - maximum number of values
133653acd3b1SBarry Smith 
133753acd3b1SBarry Smith    Output Parameter:
133853acd3b1SBarry Smith +  value - location to copy values
1339f8a50e2bSBarry Smith .  n - actual number of values found
134053acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
134153acd3b1SBarry Smith 
134253acd3b1SBarry Smith    Level: beginner
134353acd3b1SBarry Smith 
134453acd3b1SBarry Smith    Notes:
1345b32a342fSShri Abhyankar    The array can be passed as
1346b32a342fSShri Abhyankar    a comma seperated list:                                 0,1,2,3,4,5,6,7
13470fd488f5SShri Abhyankar    a range (start-end+1):                                  0-8
13480fd488f5SShri Abhyankar    a range with given increment (start-end+1:inc):         0-7:2
13490fd488f5SShri Abhyankar    a combination of values and ranges seperated by commas: 0,1-8,8-15:2
1350b32a342fSShri Abhyankar 
1351b32a342fSShri Abhyankar    There must be no intervening spaces between the values.
135253acd3b1SBarry Smith 
135353acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
135453acd3b1SBarry Smith 
1355b32a342fSShri Abhyankar    Concepts: options database^array of ints
135653acd3b1SBarry Smith 
135753acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1358acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
135953acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
136053acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1361acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1362a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsRealArray()
136353acd3b1SBarry Smith @*/
13647087cfbeSBarry Smith PetscErrorCode  PetscOptionsIntArray(const char opt[],const char text[],const char man[],PetscInt value[],PetscInt *n,PetscBool  *set)
136553acd3b1SBarry Smith {
136653acd3b1SBarry Smith   PetscErrorCode ierr;
136753acd3b1SBarry Smith   PetscInt       i;
1368e26ddf31SBarry Smith   PetscOptions   amsopt;
136953acd3b1SBarry Smith 
137053acd3b1SBarry Smith   PetscFunctionBegin;
1371b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
1372e26ddf31SBarry Smith     PetscInt *vals;
1373e26ddf31SBarry Smith 
1374e26ddf31SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_INT_ARRAY,&amsopt);CHKERRQ(ierr);
13756e02237eSPeter Brune     ierr = PetscMalloc1((*n),(PetscInt**)&amsopt->data);CHKERRQ(ierr);
1376e26ddf31SBarry Smith     vals = (PetscInt*)amsopt->data;
1377e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1378e26ddf31SBarry Smith     amsopt->arraylength = *n;
1379e26ddf31SBarry Smith   }
138053acd3b1SBarry Smith   ierr = PetscOptionsGetIntArray(PetscOptionsObject.prefix,opt,value,n,set);CHKERRQ(ierr);
138161b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
138253acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%d",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,value[0]);CHKERRQ(ierr);
138353acd3b1SBarry Smith     for (i=1; i<*n; i++) {
138453acd3b1SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,",%d",value[i]);CHKERRQ(ierr);
138553acd3b1SBarry Smith     }
13862aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
138753acd3b1SBarry Smith   }
138853acd3b1SBarry Smith   PetscFunctionReturn(0);
138953acd3b1SBarry Smith }
139053acd3b1SBarry Smith 
139153acd3b1SBarry Smith #undef __FUNCT__
139253acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsStringArray"
139353acd3b1SBarry Smith /*@C
139453acd3b1SBarry Smith    PetscOptionsStringArray - Gets an array of string values for a particular
139553acd3b1SBarry Smith    option in the database. The values must be separated with commas with
139653acd3b1SBarry Smith    no intervening spaces.
139753acd3b1SBarry Smith 
13983f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
139953acd3b1SBarry Smith 
140053acd3b1SBarry Smith    Input Parameters:
140153acd3b1SBarry Smith +  opt - the option one is seeking
140253acd3b1SBarry Smith .  text - short string describing option
140353acd3b1SBarry Smith .  man - manual page for option
140453acd3b1SBarry Smith -  nmax - maximum number of strings
140553acd3b1SBarry Smith 
140653acd3b1SBarry Smith    Output Parameter:
140753acd3b1SBarry Smith +  value - location to copy strings
140853acd3b1SBarry Smith .  nmax - actual number of strings found
140953acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
141053acd3b1SBarry Smith 
141153acd3b1SBarry Smith    Level: beginner
141253acd3b1SBarry Smith 
141353acd3b1SBarry Smith    Notes:
141453acd3b1SBarry Smith    The user should pass in an array of pointers to char, to hold all the
141553acd3b1SBarry Smith    strings returned by this function.
141653acd3b1SBarry Smith 
141753acd3b1SBarry Smith    The user is responsible for deallocating the strings that are
141853acd3b1SBarry Smith    returned. The Fortran interface for this routine is not supported.
141953acd3b1SBarry Smith 
142053acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
142153acd3b1SBarry Smith 
142253acd3b1SBarry Smith    Concepts: options database^array of strings
142353acd3b1SBarry Smith 
142453acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1425acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
142653acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
142753acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1428acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1429a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
143053acd3b1SBarry Smith @*/
14317087cfbeSBarry Smith PetscErrorCode  PetscOptionsStringArray(const char opt[],const char text[],const char man[],char *value[],PetscInt *nmax,PetscBool  *set)
143253acd3b1SBarry Smith {
143353acd3b1SBarry Smith   PetscErrorCode ierr;
14341ae3d29cSBarry Smith   PetscOptions   amsopt;
143553acd3b1SBarry Smith 
143653acd3b1SBarry Smith   PetscFunctionBegin;
14371ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
14381ae3d29cSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_STRING_ARRAY,&amsopt);CHKERRQ(ierr);
14396e02237eSPeter Brune     ierr = PetscMalloc1((*nmax),(char**)&amsopt->data);CHKERRQ(ierr);
1440a297a907SKarl Rupp 
14411ae3d29cSBarry Smith     amsopt->arraylength = *nmax;
14421ae3d29cSBarry Smith   }
144353acd3b1SBarry Smith   ierr = PetscOptionsGetStringArray(PetscOptionsObject.prefix,opt,value,nmax,set);CHKERRQ(ierr);
144461b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
14452aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <string1,string2,...>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
144653acd3b1SBarry Smith   }
144753acd3b1SBarry Smith   PetscFunctionReturn(0);
144853acd3b1SBarry Smith }
144953acd3b1SBarry Smith 
1450e2446a98SMatthew Knepley #undef __FUNCT__
1451acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolArray"
1452e2446a98SMatthew Knepley /*@C
1453acfcf0e5SJed Brown    PetscOptionsBoolArray - Gets an array of logical values (true or false) for a particular
1454e2446a98SMatthew Knepley    option in the database. The values must be separated with commas with
1455e2446a98SMatthew Knepley    no intervening spaces.
1456e2446a98SMatthew Knepley 
14573f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
1458e2446a98SMatthew Knepley 
1459e2446a98SMatthew Knepley    Input Parameters:
1460e2446a98SMatthew Knepley +  opt - the option one is seeking
1461e2446a98SMatthew Knepley .  text - short string describing option
1462e2446a98SMatthew Knepley .  man - manual page for option
1463e2446a98SMatthew Knepley -  nmax - maximum number of values
1464e2446a98SMatthew Knepley 
1465e2446a98SMatthew Knepley    Output Parameter:
1466e2446a98SMatthew Knepley +  value - location to copy values
1467e2446a98SMatthew Knepley .  nmax - actual number of values found
1468e2446a98SMatthew Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
1469e2446a98SMatthew Knepley 
1470e2446a98SMatthew Knepley    Level: beginner
1471e2446a98SMatthew Knepley 
1472e2446a98SMatthew Knepley    Notes:
1473e2446a98SMatthew Knepley    The user should pass in an array of doubles
1474e2446a98SMatthew Knepley 
1475e2446a98SMatthew Knepley    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
1476e2446a98SMatthew Knepley 
1477e2446a98SMatthew Knepley    Concepts: options database^array of strings
1478e2446a98SMatthew Knepley 
1479e2446a98SMatthew Knepley .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1480acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
1481e2446a98SMatthew Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1482e2446a98SMatthew Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1483acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1484a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
1485e2446a98SMatthew Knepley @*/
14867087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolArray(const char opt[],const char text[],const char man[],PetscBool value[],PetscInt *n,PetscBool *set)
1487e2446a98SMatthew Knepley {
1488e2446a98SMatthew Knepley   PetscErrorCode ierr;
1489e2446a98SMatthew Knepley   PetscInt       i;
14901ae3d29cSBarry Smith   PetscOptions   amsopt;
1491e2446a98SMatthew Knepley 
1492e2446a98SMatthew Knepley   PetscFunctionBegin;
14931ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
1494ace3abfcSBarry Smith     PetscBool *vals;
14951ae3d29cSBarry Smith 
14967781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL_ARRAY,&amsopt);CHKERRQ(ierr);
14976e02237eSPeter Brune     ierr = PetscMalloc1((*n),(PetscBool**)&amsopt->data);CHKERRQ(ierr);
1498ace3abfcSBarry Smith     vals = (PetscBool*)amsopt->data;
14991ae3d29cSBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
15001ae3d29cSBarry Smith     amsopt->arraylength = *n;
15011ae3d29cSBarry Smith   }
1502acfcf0e5SJed Brown   ierr = PetscOptionsGetBoolArray(PetscOptionsObject.prefix,opt,value,n,set);CHKERRQ(ierr);
1503e2446a98SMatthew Knepley   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
1504e2446a98SMatthew Knepley     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%d",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,value[0]);CHKERRQ(ierr);
1505e2446a98SMatthew Knepley     for (i=1; i<*n; i++) {
1506e2446a98SMatthew Knepley       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,",%d",value[i]);CHKERRQ(ierr);
1507e2446a98SMatthew Knepley     }
15082aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
1509e2446a98SMatthew Knepley   }
1510e2446a98SMatthew Knepley   PetscFunctionReturn(0);
1511e2446a98SMatthew Knepley }
1512e2446a98SMatthew Knepley 
15138cc676e6SMatthew G Knepley #undef __FUNCT__
15148cc676e6SMatthew G Knepley #define __FUNCT__ "PetscOptionsViewer"
15158cc676e6SMatthew G Knepley /*@C
1516d1da0b69SBarry Smith    PetscOptionsViewer - Gets a viewer appropriate for the type indicated by the user
15178cc676e6SMatthew G Knepley 
15188cc676e6SMatthew G Knepley    Logically Collective on the communicator passed in PetscOptionsBegin()
15198cc676e6SMatthew G Knepley 
15208cc676e6SMatthew G Knepley    Input Parameters:
15218cc676e6SMatthew G Knepley +  opt - option name
15228cc676e6SMatthew G Knepley .  text - short string that describes the option
15238cc676e6SMatthew G Knepley -  man - manual page with additional information on option
15248cc676e6SMatthew G Knepley 
15258cc676e6SMatthew G Knepley    Output Parameter:
15268cc676e6SMatthew G Knepley +  viewer - the viewer
15278cc676e6SMatthew G Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
15288cc676e6SMatthew G Knepley 
15298cc676e6SMatthew G Knepley    Level: beginner
15308cc676e6SMatthew G Knepley 
15318cc676e6SMatthew G Knepley    Concepts: options database^has int
15328cc676e6SMatthew G Knepley 
15338cc676e6SMatthew G Knepley    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
15348cc676e6SMatthew G Knepley 
1535d1da0b69SBarry Smith    See PetscOptionsGetVieweer() for the format of the supplied viewer and its options
15368cc676e6SMatthew G Knepley 
15378cc676e6SMatthew G Knepley .seealso: PetscOptionsGetViewer(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
15388cc676e6SMatthew G Knepley           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
15398cc676e6SMatthew G Knepley           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
15408cc676e6SMatthew G Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
15418cc676e6SMatthew G Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
15428cc676e6SMatthew G Knepley           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1543a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
15448cc676e6SMatthew G Knepley @*/
1545cffb1e40SBarry Smith PetscErrorCode  PetscOptionsViewer(const char opt[],const char text[],const char man[],PetscViewer *viewer,PetscViewerFormat *format,PetscBool  *set)
15468cc676e6SMatthew G Knepley {
15478cc676e6SMatthew G Knepley   PetscErrorCode ierr;
15488cc676e6SMatthew G Knepley   PetscOptions   amsopt;
15498cc676e6SMatthew G Knepley 
15508cc676e6SMatthew G Knepley   PetscFunctionBegin;
15518cc676e6SMatthew G Knepley   if (!PetscOptionsPublishCount) {
15528cc676e6SMatthew G Knepley     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
155364facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
15545b02f95dSBarry Smith     ierr = PetscStrdup("",(char**)&amsopt->data);CHKERRQ(ierr);
15558cc676e6SMatthew G Knepley   }
1556cffb1e40SBarry Smith   ierr = PetscOptionsGetViewer(PetscOptionsObject.comm,PetscOptionsObject.prefix,opt,viewer,format,set);CHKERRQ(ierr);
15578cc676e6SMatthew G Knepley   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
15588cc676e6SMatthew G Knepley     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%s>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,"",text,ManSection(man));CHKERRQ(ierr);
15598cc676e6SMatthew G Knepley   }
15608cc676e6SMatthew G Knepley   PetscFunctionReturn(0);
15618cc676e6SMatthew G Knepley }
15628cc676e6SMatthew G Knepley 
156353acd3b1SBarry Smith 
156453acd3b1SBarry Smith #undef __FUNCT__
156553acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsHead"
156653acd3b1SBarry Smith /*@C
1567b52f573bSBarry Smith      PetscOptionsHead - Puts a heading before listing any more published options. Used, for example,
156853acd3b1SBarry Smith             in KSPSetFromOptions_GMRES().
156953acd3b1SBarry Smith 
15703f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
157153acd3b1SBarry Smith 
157253acd3b1SBarry Smith    Input Parameter:
157353acd3b1SBarry Smith .   head - the heading text
157453acd3b1SBarry Smith 
157553acd3b1SBarry Smith 
157653acd3b1SBarry Smith    Level: intermediate
157753acd3b1SBarry Smith 
157853acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
157953acd3b1SBarry Smith 
1580b52f573bSBarry Smith           Can be followed by a call to PetscOptionsTail() in the same function.
158153acd3b1SBarry Smith 
158253acd3b1SBarry Smith    Concepts: options database^subheading
158353acd3b1SBarry Smith 
158453acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1585acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
158653acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
158753acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1588acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1589a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
159053acd3b1SBarry Smith @*/
15917087cfbeSBarry Smith PetscErrorCode  PetscOptionsHead(const char head[])
159253acd3b1SBarry Smith {
159353acd3b1SBarry Smith   PetscErrorCode ierr;
159453acd3b1SBarry Smith 
159553acd3b1SBarry Smith   PetscFunctionBegin;
159661b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
159753acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  %s\n",head);CHKERRQ(ierr);
159853acd3b1SBarry Smith   }
159953acd3b1SBarry Smith   PetscFunctionReturn(0);
160053acd3b1SBarry Smith }
160153acd3b1SBarry Smith 
160253acd3b1SBarry Smith 
160353acd3b1SBarry Smith 
160453acd3b1SBarry Smith 
160553acd3b1SBarry Smith 
160653acd3b1SBarry Smith 
1607