xref: /petsc/src/sys/objects/aoptions.c (revision a23eb9822ea42bf6e58c1b09a9c02966e33aea38)
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 
149ead66b60SBarry Smith #undef __FUNCT__
150ead66b60SBarry Smith #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 
3769c1e0ce8SBarry Smith static const char *OptionsHeader = "<head>\n"
37723a1ff79SBarry Smith                                    "<script type=\"text/javascript\" src=\"http://www.mcs.anl.gov/research/projects/saws/js/jquery-1.9.1.js\"></script>\n"
37823a1ff79SBarry Smith                                    "<script type=\"text/javascript\" src=\"http://www.mcs.anl.gov/research/projects/saws/js/SAWs.js\"></script>\n"
379d1fc0251SBarry Smith                                    "<script type=\"text/javascript\" src=\"js/PETSc.js\"></script>\n"
38064bbc9efSBarry Smith                                    "<script>\n"
38164bbc9efSBarry Smith                                       "jQuery(document).ready(function() {\n"
38264bbc9efSBarry Smith                                          "PETSc.getAndDisplayDirectory(null,\"#variablesInfo\")\n"
38364bbc9efSBarry Smith                                       "})\n"
38464bbc9efSBarry Smith                                   "</script>\n"
38564bbc9efSBarry Smith                                   "</head>\n";
386426f447eSSurtai Han static const char *OptionsBodyBottom = "<div id=\"variablesInfo\" style=\"height:550px;overflow:scroll;\"></div>\n<br>\n</body>";
38764bbc9efSBarry Smith 
388d5649816SBarry Smith #undef __FUNCT__
3897781c08eSBarry Smith #define __FUNCT__ "PetscOptionsSAWsInput"
390b3506946SBarry Smith /*
3917781c08eSBarry Smith     PetscOptionsSAWsInput - Presents all the PETSc Options processed by the program so the user may change them at runtime using the SAWs
392b3506946SBarry Smith 
393b3506946SBarry Smith     Bugs:
394b3506946SBarry Smith +    All processes must traverse through the exact same set of option queries do to the call to PetscScanString()
395b3506946SBarry Smith .    Internal strings have arbitrary length and string copies are not checked that they fit into string space
396b3506946SBarry Smith -    Only works for PetscInt == int, PetscReal == double etc
397b3506946SBarry Smith 
398b3506946SBarry Smith 
399b3506946SBarry Smith */
400475446a1SBarry Smith PetscErrorCode PetscOptionsSAWsInput()
401b3506946SBarry Smith {
402b3506946SBarry Smith   PetscErrorCode ierr;
403b3506946SBarry Smith   PetscOptions   next     = PetscOptionsObject.next;
404d5649816SBarry Smith   static int     mancount = 0;
405b3506946SBarry Smith   char           options[16];
4067aab2a10SBarry Smith   PetscBool      changedmethod = PETSC_FALSE;
407*a23eb982SSurtai Han   PetscBool      stopasking    = PETSC_FALSE;
40888a9752cSBarry Smith   char           manname[16],textname[16];
4092657e9d9SBarry Smith   char           dir[1024];
410b3506946SBarry Smith 
411b3506946SBarry Smith   /* the next line is a bug, this will only work if all processors are here, the comm passed in is ignored!!! */
412b3506946SBarry Smith   sprintf(options,"Options_%d",count++);
413a297a907SKarl Rupp 
414e04113cfSBarry Smith   PetscOptionsObject.pprefix = PetscOptionsObject.prefix; /* SAWs will change this, so cannot pass prefix directly */
4151bc75a8dSBarry Smith 
416d50831c4SBarry Smith   ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s","_title");CHKERRQ(ierr);
4177781c08eSBarry Smith   PetscStackCallSAWs(SAWs_Register,(dir,&PetscOptionsObject.title,1,SAWs_READ,SAWs_STRING));
4187781c08eSBarry Smith   ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s","prefix");CHKERRQ(ierr);
4192657e9d9SBarry Smith   PetscStackCallSAWs(SAWs_Register,(dir,&PetscOptionsObject.pprefix,1,SAWs_READ,SAWs_STRING));
4202657e9d9SBarry Smith   PetscStackCallSAWs(SAWs_Register,("/PETSc/Options/ChangedMethod",&changedmethod,1,SAWs_WRITE,SAWs_BOOLEAN));
421*a23eb982SSurtai Han   PetscStackCallSAWs(SAWs_Register,("/PETSc/Options/StopAsking",&stopasking,1,SAWs_WRITE,SAWs_BOOLEAN));
422b3506946SBarry Smith 
423b3506946SBarry Smith   while (next) {
424d50831c4SBarry Smith     sprintf(manname,"_man_%d",mancount);
4252657e9d9SBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",manname);CHKERRQ(ierr);
4267781c08eSBarry Smith     PetscStackCallSAWs(SAWs_Register,(dir,&next->man,1,SAWs_READ,SAWs_STRING));
427d50831c4SBarry Smith     sprintf(textname,"_text_%d",mancount++);
4287781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",textname);CHKERRQ(ierr);
4297781c08eSBarry Smith     PetscStackCallSAWs(SAWs_Register,(dir,&next->text,1,SAWs_READ,SAWs_STRING));
4309f32e415SBarry Smith 
431b3506946SBarry Smith     switch (next->type) {
432b3506946SBarry Smith     case OPTION_HEAD:
433b3506946SBarry Smith       break;
434b3506946SBarry Smith     case OPTION_INT_ARRAY:
4357781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4362657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_INT));
437b3506946SBarry Smith       break;
438b3506946SBarry Smith     case OPTION_REAL_ARRAY:
4397781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4402657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_DOUBLE));
441b3506946SBarry Smith       break;
442b3506946SBarry Smith     case OPTION_INT:
4437781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4442657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_INT));
445b3506946SBarry Smith       break;
446b3506946SBarry Smith     case OPTION_REAL:
4477781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4482657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_DOUBLE));
449b3506946SBarry Smith       break;
4507781c08eSBarry Smith     case OPTION_BOOL:
4517781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4522657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_BOOLEAN));
4531ae3d29cSBarry Smith       break;
4547781c08eSBarry Smith     case OPTION_BOOL_ARRAY:
4557781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4562657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_BOOLEAN));
45771f08665SBarry Smith       break;
458b3506946SBarry Smith     case OPTION_STRING:
4597781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4607781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
4611ae3d29cSBarry Smith       break;
4621ae3d29cSBarry Smith     case OPTION_STRING_ARRAY:
4637781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4642657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_STRING));
465b3506946SBarry Smith       break;
466a264d7a6SBarry Smith     case OPTION_FLIST:
467a264d7a6SBarry Smith       {
468a264d7a6SBarry Smith       PetscInt ntext;
4697781c08eSBarry Smith       ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4707781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
471a264d7a6SBarry Smith       ierr = PetscFunctionListGet(next->flist,(const char***)&next->edata,&ntext);CHKERRQ(ierr);
472a264d7a6SBarry Smith       PetscStackCallSAWs(SAWs_Set_Legal_Variable_Values,(dir,ntext,next->edata));
473a264d7a6SBarry Smith       }
474a264d7a6SBarry Smith       break;
4751ae3d29cSBarry Smith     case OPTION_ELIST:
476a264d7a6SBarry Smith       {
477a264d7a6SBarry Smith       PetscInt ntext = next->nlist;
4787781c08eSBarry Smith       ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4797781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
480ead66b60SBarry Smith       ierr = PetscMalloc1((ntext+1),(char***)&next->edata);CHKERRQ(ierr);
4811ae3d29cSBarry Smith       ierr = PetscMemcpy(next->edata,next->list,ntext*sizeof(char*));CHKERRQ(ierr);
482a264d7a6SBarry Smith       PetscStackCallSAWs(SAWs_Set_Legal_Variable_Values,(dir,ntext,next->edata));
483a264d7a6SBarry Smith       }
484a264d7a6SBarry Smith       break;
485b3506946SBarry Smith     default:
486b3506946SBarry Smith       break;
487b3506946SBarry Smith     }
488b3506946SBarry Smith     next = next->next;
489b3506946SBarry Smith   }
490b3506946SBarry Smith 
491b3506946SBarry Smith   /* wait until accessor has unlocked the memory */
49264bbc9efSBarry Smith   PetscStackCallSAWs(SAWs_Push_Header,("index.html",OptionsHeader));
49364bbc9efSBarry Smith   PetscStackCallSAWs(SAWs_Push_Body,("index.html",2,OptionsBodyBottom));
4947aab2a10SBarry Smith   ierr = PetscSAWsBlock();CHKERRQ(ierr);
49564bbc9efSBarry Smith   PetscStackCallSAWs(SAWs_Pop_Header,("index.html"));
49664bbc9efSBarry Smith   PetscStackCallSAWs(SAWs_Pop_Body,("index.html",2));
497b3506946SBarry Smith 
49888a9752cSBarry Smith   /* determine if any values have been set in GUI */
49988a9752cSBarry Smith   next = PetscOptionsObject.next;
50088a9752cSBarry Smith   while (next) {
50188a9752cSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
50288a9752cSBarry Smith     PetscStackCallSAWs(SAWs_Selected,(dir,&next->set));
50388a9752cSBarry Smith     next = next->next;
50488a9752cSBarry Smith   }
50588a9752cSBarry Smith 
506b3506946SBarry Smith   /* reset counter to -2; this updates the screen with the new options for the selected method */
507b3506946SBarry Smith   if (changedmethod) PetscOptionsPublishCount = -2;
508b3506946SBarry Smith 
509*a23eb982SSurtai Han   if (stopasking) {
510*a23eb982SSurtai Han     PetscOptionsPublish      = PETSC_FALSE;
511*a23eb982SSurtai Han     PetscOptionsPublishCount = 0;//do not ask for same thing again
512*a23eb982SSurtai Han   }
513*a23eb982SSurtai Han 
5149a492a5cSBarry Smith   PetscStackCallSAWs(SAWs_Delete,("/PETSc/Options"));
515b3506946SBarry Smith   PetscFunctionReturn(0);
516b3506946SBarry Smith }
517b3506946SBarry Smith #endif
518b3506946SBarry Smith 
5196356e834SBarry Smith #undef __FUNCT__
52053acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEnd_Private"
52153acd3b1SBarry Smith PetscErrorCode PetscOptionsEnd_Private(void)
52253acd3b1SBarry Smith {
52353acd3b1SBarry Smith   PetscErrorCode ierr;
5246356e834SBarry Smith   PetscOptions   last;
5256356e834SBarry Smith   char           option[256],value[1024],tmp[32];
526330cf3c9SBarry Smith   size_t         j;
52753acd3b1SBarry Smith 
52853acd3b1SBarry Smith   PetscFunctionBegin;
529aee2cecaSBarry Smith   if (PetscOptionsObject.next) {
530b3506946SBarry Smith     if (!PetscOptionsPublishCount) {
531a264d7a6SBarry Smith #if defined(PETSC_HAVE_SAWS)
532475446a1SBarry Smith       ierr = PetscOptionsSAWsInput();CHKERRQ(ierr);
533b3506946SBarry Smith #else
53471f08665SBarry Smith       ierr = PetscOptionsGetFromTextInput();CHKERRQ(ierr);
535b3506946SBarry Smith #endif
536aee2cecaSBarry Smith     }
537aee2cecaSBarry Smith   }
5386356e834SBarry Smith 
5396356e834SBarry Smith   /* reset counter to -2; this updates the screen with the new options for the selected method */
5406356e834SBarry Smith   if (PetscOptionsObject.changedmethod) PetscOptionsPublishCount = -2;
54161b37b28SSatish Balay   /* reset alreadyprinted flag */
54261b37b28SSatish Balay   PetscOptionsObject.alreadyprinted = PETSC_FALSE;
5433194b578SJed Brown   if (PetscOptionsObject.object) PetscOptionsObject.object->optionsprinted = PETSC_TRUE;
5440298fd71SBarry Smith   PetscOptionsObject.object = NULL;
5456356e834SBarry Smith 
5466356e834SBarry Smith   while (PetscOptionsObject.next) {
5476356e834SBarry Smith     if (PetscOptionsObject.next->set) {
5486356e834SBarry Smith       if (PetscOptionsObject.prefix) {
5496356e834SBarry Smith         ierr = PetscStrcpy(option,"-");CHKERRQ(ierr);
5506356e834SBarry Smith         ierr = PetscStrcat(option,PetscOptionsObject.prefix);CHKERRQ(ierr);
5516356e834SBarry Smith         ierr = PetscStrcat(option,PetscOptionsObject.next->option+1);CHKERRQ(ierr);
5526356e834SBarry Smith       } else {
5536356e834SBarry Smith         ierr = PetscStrcpy(option,PetscOptionsObject.next->option);CHKERRQ(ierr);
5546356e834SBarry Smith       }
5556356e834SBarry Smith 
5566356e834SBarry Smith       switch (PetscOptionsObject.next->type) {
5576356e834SBarry Smith       case OPTION_HEAD:
5586356e834SBarry Smith         break;
559e26ddf31SBarry Smith       case OPTION_INT_ARRAY:
560e26ddf31SBarry Smith         sprintf(value,"%d",(int)((PetscInt*)PetscOptionsObject.next->data)[0]);
561e26ddf31SBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
562e26ddf31SBarry Smith           sprintf(tmp,"%d",(int)((PetscInt*)PetscOptionsObject.next->data)[j]);
563e26ddf31SBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
564e26ddf31SBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
565e26ddf31SBarry Smith         }
566e26ddf31SBarry Smith         break;
5676356e834SBarry Smith       case OPTION_INT:
5687a72a596SBarry Smith         sprintf(value,"%d",(int) *(PetscInt*)PetscOptionsObject.next->data);
5696356e834SBarry Smith         break;
5706356e834SBarry Smith       case OPTION_REAL:
5717a72a596SBarry Smith         sprintf(value,"%g",(double) *(PetscReal*)PetscOptionsObject.next->data);
5726356e834SBarry Smith         break;
5736356e834SBarry Smith       case OPTION_REAL_ARRAY:
5747a72a596SBarry Smith         sprintf(value,"%g",(double)((PetscReal*)PetscOptionsObject.next->data)[0]);
5756356e834SBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
5767a72a596SBarry Smith           sprintf(tmp,"%g",(double)((PetscReal*)PetscOptionsObject.next->data)[j]);
5776356e834SBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5786356e834SBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5796356e834SBarry Smith         }
5806356e834SBarry Smith         break;
5817781c08eSBarry Smith       case OPTION_BOOL:
58271f08665SBarry Smith         sprintf(value,"%d",*(int*)PetscOptionsObject.next->data);
5836356e834SBarry Smith         break;
5847781c08eSBarry Smith       case OPTION_BOOL_ARRAY:
585ace3abfcSBarry Smith         sprintf(value,"%d",(int)((PetscBool*)PetscOptionsObject.next->data)[0]);
5861ae3d29cSBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
587ace3abfcSBarry Smith           sprintf(tmp,"%d",(int)((PetscBool*)PetscOptionsObject.next->data)[j]);
5881ae3d29cSBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5891ae3d29cSBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5901ae3d29cSBarry Smith         }
5911ae3d29cSBarry Smith         break;
592a264d7a6SBarry Smith       case OPTION_FLIST:
5931ae3d29cSBarry Smith       case OPTION_ELIST:
5947781c08eSBarry Smith         ierr = PetscStrcpy(value,(char*)PetscOptionsObject.next->data);CHKERRQ(ierr);
5956356e834SBarry Smith         break;
5961ae3d29cSBarry Smith       case OPTION_STRING:
597475446a1SBarry Smith         ierr = PetscStrcpy(value,(char*)PetscOptionsObject.next->data);CHKERRQ(ierr);
5988da2146eSBarry Smith         break;
5991ae3d29cSBarry Smith       case OPTION_STRING_ARRAY:
6001ae3d29cSBarry Smith         sprintf(value,"%s",((char**)PetscOptionsObject.next->data)[0]);
6011ae3d29cSBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
6021ae3d29cSBarry Smith           sprintf(tmp,"%s",((char**)PetscOptionsObject.next->data)[j]);
6031ae3d29cSBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
6041ae3d29cSBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
6051ae3d29cSBarry Smith         }
6066356e834SBarry Smith         break;
6076356e834SBarry Smith       }
6086356e834SBarry Smith       ierr = PetscOptionsSetValue(option,value);CHKERRQ(ierr);
6096356e834SBarry Smith     }
610503cfb0cSBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->text);CHKERRQ(ierr);
611503cfb0cSBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->option);CHKERRQ(ierr);
6126356e834SBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->man);CHKERRQ(ierr);
61371f08665SBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->edata);CHKERRQ(ierr);
614c979a496SBarry Smith 
615c979a496SBarry Smith     if ((PetscOptionsObject.next->type == OPTION_STRING) || (PetscOptionsObject.next->type == OPTION_FLIST) || (PetscOptionsObject.next->type == OPTION_ELIST)){
61664facd6cSBarry Smith       /* must use system free since SAWs may have allocated it */
617c979a496SBarry Smith       free(PetscOptionsObject.next->data);
618c979a496SBarry Smith     } else {
6197781c08eSBarry Smith       ierr   = PetscFree(PetscOptionsObject.next->data);CHKERRQ(ierr);
620c979a496SBarry Smith     }
6217781c08eSBarry Smith 
6226356e834SBarry Smith     last                    = PetscOptionsObject.next;
6236356e834SBarry Smith     PetscOptionsObject.next = PetscOptionsObject.next->next;
6246356e834SBarry Smith     ierr                    = PetscFree(last);CHKERRQ(ierr);
6256356e834SBarry Smith   }
6266356e834SBarry Smith   PetscOptionsObject.next = 0;
62776fe2943SSurtai Han 
62876fe2943SSurtai Han   ierr = PetscFree(PetscOptionsObject.title);CHKERRQ(ierr);
62976fe2943SSurtai Han   ierr = PetscFree(PetscOptionsObject.prefix);CHKERRQ(ierr);
63053acd3b1SBarry Smith   PetscFunctionReturn(0);
63153acd3b1SBarry Smith }
63253acd3b1SBarry Smith 
63353acd3b1SBarry Smith #undef __FUNCT__
63453acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEnum"
63553acd3b1SBarry Smith /*@C
63653acd3b1SBarry Smith    PetscOptionsEnum - Gets the enum value for a particular option in the database.
63753acd3b1SBarry Smith 
6383f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
63953acd3b1SBarry Smith 
64053acd3b1SBarry Smith    Input Parameters:
64153acd3b1SBarry Smith +  opt - option name
64253acd3b1SBarry Smith .  text - short string that describes the option
64353acd3b1SBarry Smith .  man - manual page with additional information on option
64453acd3b1SBarry Smith .  list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null
64553acd3b1SBarry Smith -  defaultv - the default (current) value
64653acd3b1SBarry Smith 
64753acd3b1SBarry Smith    Output Parameter:
64853acd3b1SBarry Smith +  value - the  value to return
649b32e0204SMatthew G Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
65053acd3b1SBarry Smith 
65153acd3b1SBarry Smith    Level: beginner
65253acd3b1SBarry Smith 
65353acd3b1SBarry Smith    Concepts: options database
65453acd3b1SBarry Smith 
65553acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
65653acd3b1SBarry Smith 
65753acd3b1SBarry Smith           list is usually something like PCASMTypes or some other predefined list of enum names
65853acd3b1SBarry Smith 
65953acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
660acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
661acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
66253acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
66353acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
664acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
665a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
66653acd3b1SBarry Smith @*/
6677087cfbeSBarry Smith PetscErrorCode  PetscOptionsEnum(const char opt[],const char text[],const char man[],const char *const *list,PetscEnum defaultv,PetscEnum *value,PetscBool  *set)
66853acd3b1SBarry Smith {
66953acd3b1SBarry Smith   PetscErrorCode ierr;
67053acd3b1SBarry Smith   PetscInt       ntext = 0;
671aa5bb8c0SSatish Balay   PetscInt       tval;
672ace3abfcSBarry Smith   PetscBool      tflg;
67353acd3b1SBarry Smith 
67453acd3b1SBarry Smith   PetscFunctionBegin;
67553acd3b1SBarry Smith   while (list[ntext++]) {
676e32f2f54SBarry Smith     if (ntext > 50) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument appears to be wrong or have more than 50 entries");
67753acd3b1SBarry Smith   }
678e32f2f54SBarry Smith   if (ntext < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument must have at least two entries: typename and type prefix");
67953acd3b1SBarry Smith   ntext -= 3;
680aa5bb8c0SSatish Balay   ierr   = PetscOptionsEList(opt,text,man,list,ntext,list[defaultv],&tval,&tflg);CHKERRQ(ierr);
681aa5bb8c0SSatish Balay   /* with PETSC_USE_64BIT_INDICES sizeof(PetscInt) != sizeof(PetscEnum) */
682aa5bb8c0SSatish Balay   if (tflg) *value = (PetscEnum)tval;
683aa5bb8c0SSatish Balay   if (set)  *set   = tflg;
68453acd3b1SBarry Smith   PetscFunctionReturn(0);
68553acd3b1SBarry Smith }
68653acd3b1SBarry Smith 
68753acd3b1SBarry Smith /* -------------------------------------------------------------------------------------------------------------*/
68853acd3b1SBarry Smith #undef __FUNCT__
68953acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsInt"
69053acd3b1SBarry Smith /*@C
69153acd3b1SBarry Smith    PetscOptionsInt - Gets the integer value for a particular option in the database.
69253acd3b1SBarry Smith 
6933f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
69453acd3b1SBarry Smith 
69553acd3b1SBarry Smith    Input Parameters:
69653acd3b1SBarry Smith +  opt - option name
69753acd3b1SBarry Smith .  text - short string that describes the option
69853acd3b1SBarry Smith .  man - manual page with additional information on option
699868c398cSBarry Smith -  defaultv - the default (current) value, if the user does not provide a value this is returned in value
70053acd3b1SBarry Smith 
70153acd3b1SBarry Smith    Output Parameter:
70253acd3b1SBarry Smith +  value - the integer value to return
70353acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
70453acd3b1SBarry Smith 
70553acd3b1SBarry Smith    Level: beginner
70653acd3b1SBarry Smith 
70753acd3b1SBarry Smith    Concepts: options database^has int
70853acd3b1SBarry Smith 
70953acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
71053acd3b1SBarry Smith 
71153acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
712acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
713acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
71453acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
71553acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
716acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
717a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
71853acd3b1SBarry Smith @*/
7197087cfbeSBarry Smith PetscErrorCode  PetscOptionsInt(const char opt[],const char text[],const char man[],PetscInt defaultv,PetscInt *value,PetscBool  *set)
72053acd3b1SBarry Smith {
72153acd3b1SBarry Smith   PetscErrorCode ierr;
7226356e834SBarry Smith   PetscOptions   amsopt;
72353acd3b1SBarry Smith 
72453acd3b1SBarry Smith   PetscFunctionBegin;
725b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
7263e211508SBarry Smith     PetscBool wasset;
7273e211508SBarry Smith 
7286356e834SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_INT,&amsopt);CHKERRQ(ierr);
7296356e834SBarry Smith     ierr = PetscMalloc(sizeof(PetscInt),&amsopt->data);CHKERRQ(ierr);
7306356e834SBarry Smith     *(PetscInt*)amsopt->data = defaultv;
7313e211508SBarry Smith 
7323e211508SBarry Smith     ierr = PetscOptionsGetInt(PetscOptionsObject.prefix,opt,&defaultv,&wasset);CHKERRQ(ierr);
7333e211508SBarry Smith     if (wasset) {
7343e211508SBarry Smith       *(PetscInt*)amsopt->data = defaultv;
7353e211508SBarry Smith     }
736af6d86caSBarry Smith   }
73753acd3b1SBarry Smith   ierr = PetscOptionsGetInt(PetscOptionsObject.prefix,opt,value,set);CHKERRQ(ierr);
73861b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
7392aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%d>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,defaultv,text,ManSection(man));CHKERRQ(ierr);
74053acd3b1SBarry Smith   }
74153acd3b1SBarry Smith   PetscFunctionReturn(0);
74253acd3b1SBarry Smith }
74353acd3b1SBarry Smith 
74453acd3b1SBarry Smith #undef __FUNCT__
74553acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsString"
74653acd3b1SBarry Smith /*@C
74753acd3b1SBarry Smith    PetscOptionsString - Gets the string value for a particular option in the database.
74853acd3b1SBarry Smith 
7493f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
75053acd3b1SBarry Smith 
75153acd3b1SBarry Smith    Input Parameters:
75253acd3b1SBarry Smith +  opt - option name
75353acd3b1SBarry Smith .  text - short string that describes the option
75453acd3b1SBarry Smith .  man - manual page with additional information on option
755868c398cSBarry Smith .  defaultv - the default (current) value, if the user does not provide a value this is returned in value
756bcbf2dc5SJed Brown -  len - length of the result string including null terminator
75753acd3b1SBarry Smith 
75853acd3b1SBarry Smith    Output Parameter:
75953acd3b1SBarry Smith +  value - the value to return
76053acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
76153acd3b1SBarry Smith 
76253acd3b1SBarry Smith    Level: beginner
76353acd3b1SBarry Smith 
76453acd3b1SBarry Smith    Concepts: options database^has int
76553acd3b1SBarry Smith 
76653acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
76753acd3b1SBarry Smith 
7687fccdfe4SBarry 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).
7697fccdfe4SBarry Smith 
77053acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
771acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
772acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsReal(), PetscOptionsBool(),
77353acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
77453acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
775acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
776a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
77753acd3b1SBarry Smith @*/
7787087cfbeSBarry Smith PetscErrorCode  PetscOptionsString(const char opt[],const char text[],const char man[],const char defaultv[],char value[],size_t len,PetscBool  *set)
77953acd3b1SBarry Smith {
78053acd3b1SBarry Smith   PetscErrorCode ierr;
781aee2cecaSBarry Smith   PetscOptions   amsopt;
78253acd3b1SBarry Smith 
78353acd3b1SBarry Smith   PetscFunctionBegin;
784b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
785aee2cecaSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
78664facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
7875b02f95dSBarry Smith     ierr = PetscStrdup(defaultv ? defaultv : "",(char**)&amsopt->data);CHKERRQ(ierr);
788af6d86caSBarry Smith   }
78953acd3b1SBarry Smith   ierr = PetscOptionsGetString(PetscOptionsObject.prefix,opt,value,len,set);CHKERRQ(ierr);
79061b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
7912aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%s>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,defaultv,text,ManSection(man));CHKERRQ(ierr);
79253acd3b1SBarry Smith   }
79353acd3b1SBarry Smith   PetscFunctionReturn(0);
79453acd3b1SBarry Smith }
79553acd3b1SBarry Smith 
79653acd3b1SBarry Smith #undef __FUNCT__
79753acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsReal"
79853acd3b1SBarry Smith /*@C
79953acd3b1SBarry Smith    PetscOptionsReal - Gets the PetscReal value for a particular option in the database.
80053acd3b1SBarry Smith 
8013f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
80253acd3b1SBarry Smith 
80353acd3b1SBarry Smith    Input Parameters:
80453acd3b1SBarry Smith +  opt - option name
80553acd3b1SBarry Smith .  text - short string that describes the option
80653acd3b1SBarry Smith .  man - manual page with additional information on option
807868c398cSBarry Smith -  defaultv - the default (current) value, if the user does not provide a value this is returned in value
80853acd3b1SBarry Smith 
80953acd3b1SBarry Smith    Output Parameter:
81053acd3b1SBarry Smith +  value - the value to return
81153acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
81253acd3b1SBarry Smith 
81353acd3b1SBarry Smith    Level: beginner
81453acd3b1SBarry Smith 
81553acd3b1SBarry Smith    Concepts: options database^has int
81653acd3b1SBarry Smith 
81753acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
81853acd3b1SBarry Smith 
81953acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
820acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
821acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
82253acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
82353acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
824acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
825a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
82653acd3b1SBarry Smith @*/
8277087cfbeSBarry Smith PetscErrorCode  PetscOptionsReal(const char opt[],const char text[],const char man[],PetscReal defaultv,PetscReal *value,PetscBool  *set)
82853acd3b1SBarry Smith {
82953acd3b1SBarry Smith   PetscErrorCode ierr;
830538aa990SBarry Smith   PetscOptions   amsopt;
83153acd3b1SBarry Smith 
83253acd3b1SBarry Smith   PetscFunctionBegin;
833b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
834538aa990SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_REAL,&amsopt);CHKERRQ(ierr);
835538aa990SBarry Smith     ierr = PetscMalloc(sizeof(PetscReal),&amsopt->data);CHKERRQ(ierr);
836a297a907SKarl Rupp 
837538aa990SBarry Smith     *(PetscReal*)amsopt->data = defaultv;
838538aa990SBarry Smith   }
83953acd3b1SBarry Smith   ierr = PetscOptionsGetReal(PetscOptionsObject.prefix,opt,value,set);CHKERRQ(ierr);
84061b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
84157622a8eSBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%g>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,(double)defaultv,text,ManSection(man));CHKERRQ(ierr);
84253acd3b1SBarry Smith   }
84353acd3b1SBarry Smith   PetscFunctionReturn(0);
84453acd3b1SBarry Smith }
84553acd3b1SBarry Smith 
84653acd3b1SBarry Smith #undef __FUNCT__
84753acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsScalar"
84853acd3b1SBarry Smith /*@C
84953acd3b1SBarry Smith    PetscOptionsScalar - Gets the scalar value for a particular option in the database.
85053acd3b1SBarry Smith 
8513f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
85253acd3b1SBarry Smith 
85353acd3b1SBarry Smith    Input Parameters:
85453acd3b1SBarry Smith +  opt - option name
85553acd3b1SBarry Smith .  text - short string that describes the option
85653acd3b1SBarry Smith .  man - manual page with additional information on option
857868c398cSBarry Smith -  defaultv - the default (current) value, if the user does not provide a value this is returned in value
85853acd3b1SBarry Smith 
85953acd3b1SBarry Smith    Output Parameter:
86053acd3b1SBarry Smith +  value - the value to return
86153acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
86253acd3b1SBarry Smith 
86353acd3b1SBarry Smith    Level: beginner
86453acd3b1SBarry Smith 
86553acd3b1SBarry Smith    Concepts: options database^has int
86653acd3b1SBarry Smith 
86753acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
86853acd3b1SBarry Smith 
86953acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
870acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
871acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
87253acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
87353acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
874acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
875a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
87653acd3b1SBarry Smith @*/
8777087cfbeSBarry Smith PetscErrorCode  PetscOptionsScalar(const char opt[],const char text[],const char man[],PetscScalar defaultv,PetscScalar *value,PetscBool  *set)
87853acd3b1SBarry Smith {
87953acd3b1SBarry Smith   PetscErrorCode ierr;
88053acd3b1SBarry Smith 
88153acd3b1SBarry Smith   PetscFunctionBegin;
88253acd3b1SBarry Smith #if !defined(PETSC_USE_COMPLEX)
88353acd3b1SBarry Smith   ierr = PetscOptionsReal(opt,text,man,defaultv,value,set);CHKERRQ(ierr);
88453acd3b1SBarry Smith #else
88553acd3b1SBarry Smith   ierr = PetscOptionsGetScalar(PetscOptionsObject.prefix,opt,value,set);CHKERRQ(ierr);
88653acd3b1SBarry Smith #endif
88753acd3b1SBarry Smith   PetscFunctionReturn(0);
88853acd3b1SBarry Smith }
88953acd3b1SBarry Smith 
89053acd3b1SBarry Smith #undef __FUNCT__
89153acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsName"
89253acd3b1SBarry Smith /*@C
89390d69ab7SBarry 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
89490d69ab7SBarry Smith                       its value is set to false.
89553acd3b1SBarry Smith 
8963f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
89753acd3b1SBarry Smith 
89853acd3b1SBarry Smith    Input Parameters:
89953acd3b1SBarry Smith +  opt - option name
90053acd3b1SBarry Smith .  text - short string that describes the option
90153acd3b1SBarry Smith -  man - manual page with additional information on option
90253acd3b1SBarry Smith 
90353acd3b1SBarry Smith    Output Parameter:
90453acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
90553acd3b1SBarry Smith 
90653acd3b1SBarry Smith    Level: beginner
90753acd3b1SBarry Smith 
90853acd3b1SBarry Smith    Concepts: options database^has int
90953acd3b1SBarry Smith 
91053acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
91153acd3b1SBarry Smith 
91253acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
913acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
914acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
91553acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
91653acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
917acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
918a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
91953acd3b1SBarry Smith @*/
9207087cfbeSBarry Smith PetscErrorCode  PetscOptionsName(const char opt[],const char text[],const char man[],PetscBool  *flg)
92153acd3b1SBarry Smith {
92253acd3b1SBarry Smith   PetscErrorCode ierr;
9231ae3d29cSBarry Smith   PetscOptions   amsopt;
92453acd3b1SBarry Smith 
92553acd3b1SBarry Smith   PetscFunctionBegin;
9261ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
9277781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
928ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
929a297a907SKarl Rupp 
930ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
9311ae3d29cSBarry Smith   }
93253acd3b1SBarry Smith   ierr = PetscOptionsHasName(PetscOptionsObject.prefix,opt,flg);CHKERRQ(ierr);
93361b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
9342aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
93553acd3b1SBarry Smith   }
93653acd3b1SBarry Smith   PetscFunctionReturn(0);
93753acd3b1SBarry Smith }
93853acd3b1SBarry Smith 
93953acd3b1SBarry Smith #undef __FUNCT__
940a264d7a6SBarry Smith #define __FUNCT__ "PetscOptionsFList"
94153acd3b1SBarry Smith /*@C
942a264d7a6SBarry Smith      PetscOptionsFList - Puts a list of option values that a single one may be selected from
94353acd3b1SBarry Smith 
9443f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
94553acd3b1SBarry Smith 
94653acd3b1SBarry Smith    Input Parameters:
94753acd3b1SBarry Smith +  opt - option name
94853acd3b1SBarry Smith .  text - short string that describes the option
94953acd3b1SBarry Smith .  man - manual page with additional information on option
95053acd3b1SBarry Smith .  list - the possible choices
951868c398cSBarry Smith .  defaultv - the default (current) value, if the user does not provide a value this is returned in value
9523cc1e11dSBarry Smith -  len - the length of the character array value
95353acd3b1SBarry Smith 
95453acd3b1SBarry Smith    Output Parameter:
95553acd3b1SBarry Smith +  value - the value to return
95653acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
95753acd3b1SBarry Smith 
95853acd3b1SBarry Smith    Level: intermediate
95953acd3b1SBarry Smith 
96053acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
96153acd3b1SBarry Smith 
96253acd3b1SBarry Smith    See PetscOptionsEList() for when the choices are given in a string array
96353acd3b1SBarry Smith 
96453acd3b1SBarry Smith    To get a listing of all currently specified options,
96588c29154SBarry Smith     see PetscOptionsView() or PetscOptionsGetAll()
96653acd3b1SBarry Smith 
967eabe10d7SBarry Smith    Developer Note: This cannot check for invalid selection because of things like MATAIJ that are not included in the list
968eabe10d7SBarry Smith 
96953acd3b1SBarry Smith    Concepts: options database^list
97053acd3b1SBarry Smith 
97153acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
972acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
97353acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
97453acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
975acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
976a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsEnum()
97753acd3b1SBarry Smith @*/
978a264d7a6SBarry Smith PetscErrorCode  PetscOptionsFList(const char opt[],const char ltext[],const char man[],PetscFunctionList list,const char defaultv[],char value[],size_t len,PetscBool  *set)
97953acd3b1SBarry Smith {
98053acd3b1SBarry Smith   PetscErrorCode ierr;
9813cc1e11dSBarry Smith   PetscOptions   amsopt;
98253acd3b1SBarry Smith 
98353acd3b1SBarry Smith   PetscFunctionBegin;
984b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
985a264d7a6SBarry Smith     ierr = PetscOptionsCreate_Private(opt,ltext,man,OPTION_FLIST,&amsopt);CHKERRQ(ierr);
98664facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
9875b02f95dSBarry Smith     ierr = PetscStrdup(defaultv ? defaultv : "",(char**)&amsopt->data);CHKERRQ(ierr);
9883cc1e11dSBarry Smith     amsopt->flist = list;
9893cc1e11dSBarry Smith   }
99053acd3b1SBarry Smith   ierr = PetscOptionsGetString(PetscOptionsObject.prefix,opt,value,len,set);CHKERRQ(ierr);
99161b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
992140e18c1SBarry Smith     ierr = PetscFunctionListPrintTypes(PetscOptionsObject.comm,stdout,PetscOptionsObject.prefix,opt,ltext,man,list,defaultv);CHKERRQ(ierr);CHKERRQ(ierr);
99353acd3b1SBarry Smith   }
99453acd3b1SBarry Smith   PetscFunctionReturn(0);
99553acd3b1SBarry Smith }
99653acd3b1SBarry Smith 
99753acd3b1SBarry Smith #undef __FUNCT__
99853acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEList"
99953acd3b1SBarry Smith /*@C
100053acd3b1SBarry Smith      PetscOptionsEList - Puts a list of option values that a single one may be selected from
100153acd3b1SBarry Smith 
10023f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
100353acd3b1SBarry Smith 
100453acd3b1SBarry Smith    Input Parameters:
100553acd3b1SBarry Smith +  opt - option name
100653acd3b1SBarry Smith .  ltext - short string that describes the option
100753acd3b1SBarry Smith .  man - manual page with additional information on option
1008a264d7a6SBarry Smith .  list - the possible choices (one of these must be selected, anything else is invalid)
100953acd3b1SBarry Smith .  ntext - number of choices
1010868c398cSBarry Smith -  defaultv - the default (current) value, if the user does not provide a value the index of defaultv is returned
101153acd3b1SBarry Smith 
101253acd3b1SBarry Smith    Output Parameter:
101353acd3b1SBarry Smith +  value - the index of the value to return
101453acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
101553acd3b1SBarry Smith 
101653acd3b1SBarry Smith    Level: intermediate
101753acd3b1SBarry Smith 
101853acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
101953acd3b1SBarry Smith 
1020a264d7a6SBarry Smith    See PetscOptionsFList() for when the choices are given in a PetscFunctionList()
102153acd3b1SBarry Smith 
102253acd3b1SBarry Smith    Concepts: options database^list
102353acd3b1SBarry Smith 
102453acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1025acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
102653acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
102753acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1028acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1029a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEnum()
103053acd3b1SBarry Smith @*/
10317087cfbeSBarry Smith PetscErrorCode  PetscOptionsEList(const char opt[],const char ltext[],const char man[],const char *const *list,PetscInt ntext,const char defaultv[],PetscInt *value,PetscBool  *set)
103253acd3b1SBarry Smith {
103353acd3b1SBarry Smith   PetscErrorCode ierr;
103453acd3b1SBarry Smith   PetscInt       i;
10351ae3d29cSBarry Smith   PetscOptions   amsopt;
103653acd3b1SBarry Smith 
103753acd3b1SBarry Smith   PetscFunctionBegin;
10381ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
1039d5649816SBarry Smith     ierr = PetscOptionsCreate_Private(opt,ltext,man,OPTION_ELIST,&amsopt);CHKERRQ(ierr);
104064facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
10415b02f95dSBarry Smith     ierr = PetscStrdup(defaultv ? defaultv : "",(char**)&amsopt->data);CHKERRQ(ierr);
10421ae3d29cSBarry Smith     amsopt->list  = list;
10431ae3d29cSBarry Smith     amsopt->nlist = ntext;
10441ae3d29cSBarry Smith   }
104553acd3b1SBarry Smith   ierr = PetscOptionsGetEList(PetscOptionsObject.prefix,opt,list,ntext,value,set);CHKERRQ(ierr);
104661b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
104753acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%s> (choose one of)",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,defaultv);CHKERRQ(ierr);
104853acd3b1SBarry Smith     for (i=0; i<ntext; i++) {
104953acd3b1SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm," %s",list[i]);CHKERRQ(ierr);
105053acd3b1SBarry Smith     }
10512aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm," (%s)\n",ManSection(man));CHKERRQ(ierr);
105253acd3b1SBarry Smith   }
105353acd3b1SBarry Smith   PetscFunctionReturn(0);
105453acd3b1SBarry Smith }
105553acd3b1SBarry Smith 
105653acd3b1SBarry Smith #undef __FUNCT__
1057acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolGroupBegin"
105853acd3b1SBarry Smith /*@C
1059acfcf0e5SJed Brown      PetscOptionsBoolGroupBegin - First in a series of logical queries on the options database for
1060d5649816SBarry Smith        which at most a single value can be true.
106153acd3b1SBarry Smith 
10623f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
106353acd3b1SBarry Smith 
106453acd3b1SBarry Smith    Input Parameters:
106553acd3b1SBarry Smith +  opt - option name
106653acd3b1SBarry Smith .  text - short string that describes the option
106753acd3b1SBarry Smith -  man - manual page with additional information on option
106853acd3b1SBarry Smith 
106953acd3b1SBarry Smith    Output Parameter:
107053acd3b1SBarry Smith .  flg - whether that option was set or not
107153acd3b1SBarry Smith 
107253acd3b1SBarry Smith    Level: intermediate
107353acd3b1SBarry Smith 
107453acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
107553acd3b1SBarry Smith 
1076acfcf0e5SJed Brown    Must be followed by 0 or more PetscOptionsBoolGroup()s and PetscOptionsBoolGroupEnd()
107753acd3b1SBarry Smith 
107853acd3b1SBarry Smith     Concepts: options database^logical group
107953acd3b1SBarry Smith 
108053acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1081acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
108253acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
108353acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1084acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1085a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
108653acd3b1SBarry Smith @*/
10877087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolGroupBegin(const char opt[],const char text[],const char man[],PetscBool  *flg)
108853acd3b1SBarry Smith {
108953acd3b1SBarry Smith   PetscErrorCode ierr;
10901ae3d29cSBarry Smith   PetscOptions   amsopt;
109153acd3b1SBarry Smith 
109253acd3b1SBarry Smith   PetscFunctionBegin;
10931ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
10947781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1095ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1096a297a907SKarl Rupp 
1097ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
10981ae3d29cSBarry Smith   }
109968b16fdaSBarry Smith   *flg = PETSC_FALSE;
11000298fd71SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,NULL);CHKERRQ(ierr);
110161b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
110253acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  Pick at most one of -------------\n");CHKERRQ(ierr);
11032aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"    -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
110453acd3b1SBarry Smith   }
110553acd3b1SBarry Smith   PetscFunctionReturn(0);
110653acd3b1SBarry Smith }
110753acd3b1SBarry Smith 
110853acd3b1SBarry Smith #undef __FUNCT__
1109acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolGroup"
111053acd3b1SBarry Smith /*@C
1111acfcf0e5SJed Brown      PetscOptionsBoolGroup - One in a series of logical queries on the options database for
1112d5649816SBarry Smith        which at most a single value can be true.
111353acd3b1SBarry Smith 
11143f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
111553acd3b1SBarry Smith 
111653acd3b1SBarry Smith    Input Parameters:
111753acd3b1SBarry Smith +  opt - option name
111853acd3b1SBarry Smith .  text - short string that describes the option
111953acd3b1SBarry Smith -  man - manual page with additional information on option
112053acd3b1SBarry Smith 
112153acd3b1SBarry Smith    Output Parameter:
112253acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
112353acd3b1SBarry Smith 
112453acd3b1SBarry Smith    Level: intermediate
112553acd3b1SBarry Smith 
112653acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
112753acd3b1SBarry Smith 
1128acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin() and preceded a PetscOptionsBoolGroupEnd()
112953acd3b1SBarry Smith 
113053acd3b1SBarry Smith     Concepts: options database^logical group
113153acd3b1SBarry Smith 
113253acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1133acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
113453acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
113553acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1136acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1137a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
113853acd3b1SBarry Smith @*/
11397087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolGroup(const char opt[],const char text[],const char man[],PetscBool  *flg)
114053acd3b1SBarry Smith {
114153acd3b1SBarry Smith   PetscErrorCode ierr;
11421ae3d29cSBarry Smith   PetscOptions   amsopt;
114353acd3b1SBarry Smith 
114453acd3b1SBarry Smith   PetscFunctionBegin;
11451ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
11467781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1147ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1148a297a907SKarl Rupp 
1149ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
11501ae3d29cSBarry Smith   }
115117326d04SJed Brown   *flg = PETSC_FALSE;
11520298fd71SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,NULL);CHKERRQ(ierr);
115361b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
11542aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"    -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
115553acd3b1SBarry Smith   }
115653acd3b1SBarry Smith   PetscFunctionReturn(0);
115753acd3b1SBarry Smith }
115853acd3b1SBarry Smith 
115953acd3b1SBarry Smith #undef __FUNCT__
1160acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolGroupEnd"
116153acd3b1SBarry Smith /*@C
1162acfcf0e5SJed Brown      PetscOptionsBoolGroupEnd - Last in a series of logical queries on the options database for
1163d5649816SBarry Smith        which at most a single value can be true.
116453acd3b1SBarry Smith 
11653f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
116653acd3b1SBarry Smith 
116753acd3b1SBarry Smith    Input Parameters:
116853acd3b1SBarry Smith +  opt - option name
116953acd3b1SBarry Smith .  text - short string that describes the option
117053acd3b1SBarry Smith -  man - manual page with additional information on option
117153acd3b1SBarry Smith 
117253acd3b1SBarry Smith    Output Parameter:
117353acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
117453acd3b1SBarry Smith 
117553acd3b1SBarry Smith    Level: intermediate
117653acd3b1SBarry Smith 
117753acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
117853acd3b1SBarry Smith 
1179acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin()
118053acd3b1SBarry Smith 
118153acd3b1SBarry Smith     Concepts: options database^logical group
118253acd3b1SBarry Smith 
118353acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1184acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
118553acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
118653acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1187acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1188a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
118953acd3b1SBarry Smith @*/
11907087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolGroupEnd(const char opt[],const char text[],const char man[],PetscBool  *flg)
119153acd3b1SBarry Smith {
119253acd3b1SBarry Smith   PetscErrorCode ierr;
11931ae3d29cSBarry Smith   PetscOptions   amsopt;
119453acd3b1SBarry Smith 
119553acd3b1SBarry Smith   PetscFunctionBegin;
11961ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
11977781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1198ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1199a297a907SKarl Rupp 
1200ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
12011ae3d29cSBarry Smith   }
120217326d04SJed Brown   *flg = PETSC_FALSE;
12030298fd71SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,NULL);CHKERRQ(ierr);
120461b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
12052aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"    -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
120653acd3b1SBarry Smith   }
120753acd3b1SBarry Smith   PetscFunctionReturn(0);
120853acd3b1SBarry Smith }
120953acd3b1SBarry Smith 
121053acd3b1SBarry Smith #undef __FUNCT__
1211acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBool"
121253acd3b1SBarry Smith /*@C
1213acfcf0e5SJed Brown    PetscOptionsBool - Determines if a particular option is in the database with a true or false
121453acd3b1SBarry Smith 
12153f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
121653acd3b1SBarry Smith 
121753acd3b1SBarry Smith    Input Parameters:
121853acd3b1SBarry Smith +  opt - option name
121953acd3b1SBarry Smith .  text - short string that describes the option
1220868c398cSBarry Smith .  man - manual page with additional information on option
1221868c398cSBarry Smith -  deflt - the default value, if the user does not set a value then this value is returned in flg
122253acd3b1SBarry Smith 
122353acd3b1SBarry Smith    Output Parameter:
122453acd3b1SBarry Smith .  flg - PETSC_TRUE or PETSC_FALSE
122553acd3b1SBarry Smith .  set - PETSC_TRUE if found, else PETSC_FALSE
122653acd3b1SBarry Smith 
122753acd3b1SBarry Smith    Level: beginner
122853acd3b1SBarry Smith 
122953acd3b1SBarry Smith    Concepts: options database^logical
123053acd3b1SBarry Smith 
123153acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
123253acd3b1SBarry Smith 
123353acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
1234acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
1235acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
123653acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
123753acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1238acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1239a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
124053acd3b1SBarry Smith @*/
12417087cfbeSBarry Smith PetscErrorCode  PetscOptionsBool(const char opt[],const char text[],const char man[],PetscBool deflt,PetscBool  *flg,PetscBool  *set)
124253acd3b1SBarry Smith {
124353acd3b1SBarry Smith   PetscErrorCode ierr;
1244ace3abfcSBarry Smith   PetscBool      iset;
1245aee2cecaSBarry Smith   PetscOptions   amsopt;
124653acd3b1SBarry Smith 
124753acd3b1SBarry Smith   PetscFunctionBegin;
1248b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
12497781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1250ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1251a297a907SKarl Rupp 
1252ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = deflt;
1253af6d86caSBarry Smith   }
1254acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,&iset);CHKERRQ(ierr);
125553acd3b1SBarry Smith   if (!iset) {
125653acd3b1SBarry Smith     if (flg) *flg = deflt;
125753acd3b1SBarry Smith   }
125853acd3b1SBarry Smith   if (set) *set = iset;
125961b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
1260ace3abfcSBarry Smith     const char *v = PetscBools[deflt];
12612aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s: <%s> %s (%s)\n",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,v,text,ManSection(man));CHKERRQ(ierr);
126253acd3b1SBarry Smith   }
126353acd3b1SBarry Smith   PetscFunctionReturn(0);
126453acd3b1SBarry Smith }
126553acd3b1SBarry Smith 
126653acd3b1SBarry Smith #undef __FUNCT__
126753acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsRealArray"
126853acd3b1SBarry Smith /*@C
126953acd3b1SBarry Smith    PetscOptionsRealArray - Gets an array of double values for a particular
127053acd3b1SBarry Smith    option in the database. The values must be separated with commas with
127153acd3b1SBarry Smith    no intervening spaces.
127253acd3b1SBarry Smith 
12733f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
127453acd3b1SBarry Smith 
127553acd3b1SBarry Smith    Input Parameters:
127653acd3b1SBarry Smith +  opt - the option one is seeking
127753acd3b1SBarry Smith .  text - short string describing option
127853acd3b1SBarry Smith .  man - manual page for option
127953acd3b1SBarry Smith -  nmax - maximum number of values
128053acd3b1SBarry Smith 
128153acd3b1SBarry Smith    Output Parameter:
128253acd3b1SBarry Smith +  value - location to copy values
128353acd3b1SBarry Smith .  nmax - actual number of values found
128453acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
128553acd3b1SBarry Smith 
128653acd3b1SBarry Smith    Level: beginner
128753acd3b1SBarry Smith 
128853acd3b1SBarry Smith    Notes:
128953acd3b1SBarry Smith    The user should pass in an array of doubles
129053acd3b1SBarry Smith 
129153acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
129253acd3b1SBarry Smith 
129353acd3b1SBarry Smith    Concepts: options database^array of strings
129453acd3b1SBarry Smith 
129553acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1296acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
129753acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
129853acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1299acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1300a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
130153acd3b1SBarry Smith @*/
13027087cfbeSBarry Smith PetscErrorCode  PetscOptionsRealArray(const char opt[],const char text[],const char man[],PetscReal value[],PetscInt *n,PetscBool  *set)
130353acd3b1SBarry Smith {
130453acd3b1SBarry Smith   PetscErrorCode ierr;
130553acd3b1SBarry Smith   PetscInt       i;
1306e26ddf31SBarry Smith   PetscOptions   amsopt;
130753acd3b1SBarry Smith 
130853acd3b1SBarry Smith   PetscFunctionBegin;
1309b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
1310e26ddf31SBarry Smith     PetscReal *vals;
1311e26ddf31SBarry Smith 
1312e26ddf31SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_REAL_ARRAY,&amsopt);CHKERRQ(ierr);
13136e02237eSPeter Brune     ierr = PetscMalloc1((*n),(PetscReal**)&amsopt->data);CHKERRQ(ierr);
1314e26ddf31SBarry Smith     vals = (PetscReal*)amsopt->data;
1315e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1316e26ddf31SBarry Smith     amsopt->arraylength = *n;
1317e26ddf31SBarry Smith   }
131853acd3b1SBarry Smith   ierr = PetscOptionsGetRealArray(PetscOptionsObject.prefix,opt,value,n,set);CHKERRQ(ierr);
131961b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
132057622a8eSBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%g",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,(double)value[0]);CHKERRQ(ierr);
132153acd3b1SBarry Smith     for (i=1; i<*n; i++) {
132257622a8eSBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,",%g",(double)value[i]);CHKERRQ(ierr);
132353acd3b1SBarry Smith     }
13242aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
132553acd3b1SBarry Smith   }
132653acd3b1SBarry Smith   PetscFunctionReturn(0);
132753acd3b1SBarry Smith }
132853acd3b1SBarry Smith 
132953acd3b1SBarry Smith 
133053acd3b1SBarry Smith #undef __FUNCT__
133153acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsIntArray"
133253acd3b1SBarry Smith /*@C
133353acd3b1SBarry Smith    PetscOptionsIntArray - Gets an array of integers for a particular
1334b32a342fSShri Abhyankar    option in the database.
133553acd3b1SBarry Smith 
13363f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
133753acd3b1SBarry Smith 
133853acd3b1SBarry Smith    Input Parameters:
133953acd3b1SBarry Smith +  opt - the option one is seeking
134053acd3b1SBarry Smith .  text - short string describing option
134153acd3b1SBarry Smith .  man - manual page for option
1342f8a50e2bSBarry Smith -  n - maximum number of values
134353acd3b1SBarry Smith 
134453acd3b1SBarry Smith    Output Parameter:
134553acd3b1SBarry Smith +  value - location to copy values
1346f8a50e2bSBarry Smith .  n - actual number of values found
134753acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
134853acd3b1SBarry Smith 
134953acd3b1SBarry Smith    Level: beginner
135053acd3b1SBarry Smith 
135153acd3b1SBarry Smith    Notes:
1352b32a342fSShri Abhyankar    The array can be passed as
1353b32a342fSShri Abhyankar    a comma seperated list:                                 0,1,2,3,4,5,6,7
13540fd488f5SShri Abhyankar    a range (start-end+1):                                  0-8
13550fd488f5SShri Abhyankar    a range with given increment (start-end+1:inc):         0-7:2
13560fd488f5SShri Abhyankar    a combination of values and ranges seperated by commas: 0,1-8,8-15:2
1357b32a342fSShri Abhyankar 
1358b32a342fSShri Abhyankar    There must be no intervening spaces between the values.
135953acd3b1SBarry Smith 
136053acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
136153acd3b1SBarry Smith 
1362b32a342fSShri Abhyankar    Concepts: options database^array of ints
136353acd3b1SBarry Smith 
136453acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1365acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
136653acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
136753acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1368acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1369a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsRealArray()
137053acd3b1SBarry Smith @*/
13717087cfbeSBarry Smith PetscErrorCode  PetscOptionsIntArray(const char opt[],const char text[],const char man[],PetscInt value[],PetscInt *n,PetscBool  *set)
137253acd3b1SBarry Smith {
137353acd3b1SBarry Smith   PetscErrorCode ierr;
137453acd3b1SBarry Smith   PetscInt       i;
1375e26ddf31SBarry Smith   PetscOptions   amsopt;
137653acd3b1SBarry Smith 
137753acd3b1SBarry Smith   PetscFunctionBegin;
1378b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
1379e26ddf31SBarry Smith     PetscInt *vals;
1380e26ddf31SBarry Smith 
1381e26ddf31SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_INT_ARRAY,&amsopt);CHKERRQ(ierr);
13826e02237eSPeter Brune     ierr = PetscMalloc1((*n),(PetscInt**)&amsopt->data);CHKERRQ(ierr);
1383e26ddf31SBarry Smith     vals = (PetscInt*)amsopt->data;
1384e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1385e26ddf31SBarry Smith     amsopt->arraylength = *n;
1386e26ddf31SBarry Smith   }
138753acd3b1SBarry Smith   ierr = PetscOptionsGetIntArray(PetscOptionsObject.prefix,opt,value,n,set);CHKERRQ(ierr);
138861b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
138953acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%d",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,value[0]);CHKERRQ(ierr);
139053acd3b1SBarry Smith     for (i=1; i<*n; i++) {
139153acd3b1SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,",%d",value[i]);CHKERRQ(ierr);
139253acd3b1SBarry Smith     }
13932aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
139453acd3b1SBarry Smith   }
139553acd3b1SBarry Smith   PetscFunctionReturn(0);
139653acd3b1SBarry Smith }
139753acd3b1SBarry Smith 
139853acd3b1SBarry Smith #undef __FUNCT__
139953acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsStringArray"
140053acd3b1SBarry Smith /*@C
140153acd3b1SBarry Smith    PetscOptionsStringArray - Gets an array of string values for a particular
140253acd3b1SBarry Smith    option in the database. The values must be separated with commas with
140353acd3b1SBarry Smith    no intervening spaces.
140453acd3b1SBarry Smith 
14053f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
140653acd3b1SBarry Smith 
140753acd3b1SBarry Smith    Input Parameters:
140853acd3b1SBarry Smith +  opt - the option one is seeking
140953acd3b1SBarry Smith .  text - short string describing option
141053acd3b1SBarry Smith .  man - manual page for option
141153acd3b1SBarry Smith -  nmax - maximum number of strings
141253acd3b1SBarry Smith 
141353acd3b1SBarry Smith    Output Parameter:
141453acd3b1SBarry Smith +  value - location to copy strings
141553acd3b1SBarry Smith .  nmax - actual number of strings found
141653acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
141753acd3b1SBarry Smith 
141853acd3b1SBarry Smith    Level: beginner
141953acd3b1SBarry Smith 
142053acd3b1SBarry Smith    Notes:
142153acd3b1SBarry Smith    The user should pass in an array of pointers to char, to hold all the
142253acd3b1SBarry Smith    strings returned by this function.
142353acd3b1SBarry Smith 
142453acd3b1SBarry Smith    The user is responsible for deallocating the strings that are
142553acd3b1SBarry Smith    returned. The Fortran interface for this routine is not supported.
142653acd3b1SBarry Smith 
142753acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
142853acd3b1SBarry Smith 
142953acd3b1SBarry Smith    Concepts: options database^array of strings
143053acd3b1SBarry Smith 
143153acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1432acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
143353acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
143453acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1435acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1436a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
143753acd3b1SBarry Smith @*/
14387087cfbeSBarry Smith PetscErrorCode  PetscOptionsStringArray(const char opt[],const char text[],const char man[],char *value[],PetscInt *nmax,PetscBool  *set)
143953acd3b1SBarry Smith {
144053acd3b1SBarry Smith   PetscErrorCode ierr;
14411ae3d29cSBarry Smith   PetscOptions   amsopt;
144253acd3b1SBarry Smith 
144353acd3b1SBarry Smith   PetscFunctionBegin;
14441ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
14451ae3d29cSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_STRING_ARRAY,&amsopt);CHKERRQ(ierr);
14466e02237eSPeter Brune     ierr = PetscMalloc1((*nmax),(char**)&amsopt->data);CHKERRQ(ierr);
1447a297a907SKarl Rupp 
14481ae3d29cSBarry Smith     amsopt->arraylength = *nmax;
14491ae3d29cSBarry Smith   }
145053acd3b1SBarry Smith   ierr = PetscOptionsGetStringArray(PetscOptionsObject.prefix,opt,value,nmax,set);CHKERRQ(ierr);
145161b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
14522aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <string1,string2,...>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
145353acd3b1SBarry Smith   }
145453acd3b1SBarry Smith   PetscFunctionReturn(0);
145553acd3b1SBarry Smith }
145653acd3b1SBarry Smith 
1457e2446a98SMatthew Knepley #undef __FUNCT__
1458acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolArray"
1459e2446a98SMatthew Knepley /*@C
1460acfcf0e5SJed Brown    PetscOptionsBoolArray - Gets an array of logical values (true or false) for a particular
1461e2446a98SMatthew Knepley    option in the database. The values must be separated with commas with
1462e2446a98SMatthew Knepley    no intervening spaces.
1463e2446a98SMatthew Knepley 
14643f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
1465e2446a98SMatthew Knepley 
1466e2446a98SMatthew Knepley    Input Parameters:
1467e2446a98SMatthew Knepley +  opt - the option one is seeking
1468e2446a98SMatthew Knepley .  text - short string describing option
1469e2446a98SMatthew Knepley .  man - manual page for option
1470e2446a98SMatthew Knepley -  nmax - maximum number of values
1471e2446a98SMatthew Knepley 
1472e2446a98SMatthew Knepley    Output Parameter:
1473e2446a98SMatthew Knepley +  value - location to copy values
1474e2446a98SMatthew Knepley .  nmax - actual number of values found
1475e2446a98SMatthew Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
1476e2446a98SMatthew Knepley 
1477e2446a98SMatthew Knepley    Level: beginner
1478e2446a98SMatthew Knepley 
1479e2446a98SMatthew Knepley    Notes:
1480e2446a98SMatthew Knepley    The user should pass in an array of doubles
1481e2446a98SMatthew Knepley 
1482e2446a98SMatthew Knepley    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
1483e2446a98SMatthew Knepley 
1484e2446a98SMatthew Knepley    Concepts: options database^array of strings
1485e2446a98SMatthew Knepley 
1486e2446a98SMatthew Knepley .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1487acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
1488e2446a98SMatthew Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1489e2446a98SMatthew Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1490acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1491a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
1492e2446a98SMatthew Knepley @*/
14937087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolArray(const char opt[],const char text[],const char man[],PetscBool value[],PetscInt *n,PetscBool *set)
1494e2446a98SMatthew Knepley {
1495e2446a98SMatthew Knepley   PetscErrorCode ierr;
1496e2446a98SMatthew Knepley   PetscInt       i;
14971ae3d29cSBarry Smith   PetscOptions   amsopt;
1498e2446a98SMatthew Knepley 
1499e2446a98SMatthew Knepley   PetscFunctionBegin;
15001ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
1501ace3abfcSBarry Smith     PetscBool *vals;
15021ae3d29cSBarry Smith 
15037781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL_ARRAY,&amsopt);CHKERRQ(ierr);
15046e02237eSPeter Brune     ierr = PetscMalloc1((*n),(PetscBool**)&amsopt->data);CHKERRQ(ierr);
1505ace3abfcSBarry Smith     vals = (PetscBool*)amsopt->data;
15061ae3d29cSBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
15071ae3d29cSBarry Smith     amsopt->arraylength = *n;
15081ae3d29cSBarry Smith   }
1509acfcf0e5SJed Brown   ierr = PetscOptionsGetBoolArray(PetscOptionsObject.prefix,opt,value,n,set);CHKERRQ(ierr);
1510e2446a98SMatthew Knepley   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
1511e2446a98SMatthew Knepley     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%d",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,value[0]);CHKERRQ(ierr);
1512e2446a98SMatthew Knepley     for (i=1; i<*n; i++) {
1513e2446a98SMatthew Knepley       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,",%d",value[i]);CHKERRQ(ierr);
1514e2446a98SMatthew Knepley     }
15152aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
1516e2446a98SMatthew Knepley   }
1517e2446a98SMatthew Knepley   PetscFunctionReturn(0);
1518e2446a98SMatthew Knepley }
1519e2446a98SMatthew Knepley 
15208cc676e6SMatthew G Knepley #undef __FUNCT__
15218cc676e6SMatthew G Knepley #define __FUNCT__ "PetscOptionsViewer"
15228cc676e6SMatthew G Knepley /*@C
1523d1da0b69SBarry Smith    PetscOptionsViewer - Gets a viewer appropriate for the type indicated by the user
15248cc676e6SMatthew G Knepley 
15258cc676e6SMatthew G Knepley    Logically Collective on the communicator passed in PetscOptionsBegin()
15268cc676e6SMatthew G Knepley 
15278cc676e6SMatthew G Knepley    Input Parameters:
15288cc676e6SMatthew G Knepley +  opt - option name
15298cc676e6SMatthew G Knepley .  text - short string that describes the option
15308cc676e6SMatthew G Knepley -  man - manual page with additional information on option
15318cc676e6SMatthew G Knepley 
15328cc676e6SMatthew G Knepley    Output Parameter:
15338cc676e6SMatthew G Knepley +  viewer - the viewer
15348cc676e6SMatthew G Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
15358cc676e6SMatthew G Knepley 
15368cc676e6SMatthew G Knepley    Level: beginner
15378cc676e6SMatthew G Knepley 
15388cc676e6SMatthew G Knepley    Concepts: options database^has int
15398cc676e6SMatthew G Knepley 
15408cc676e6SMatthew G Knepley    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
15418cc676e6SMatthew G Knepley 
1542d1da0b69SBarry Smith    See PetscOptionsGetVieweer() for the format of the supplied viewer and its options
15438cc676e6SMatthew G Knepley 
15448cc676e6SMatthew G Knepley .seealso: PetscOptionsGetViewer(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
15458cc676e6SMatthew G Knepley           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
15468cc676e6SMatthew G Knepley           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
15478cc676e6SMatthew G Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
15488cc676e6SMatthew G Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
15498cc676e6SMatthew G Knepley           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1550a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
15518cc676e6SMatthew G Knepley @*/
1552cffb1e40SBarry Smith PetscErrorCode  PetscOptionsViewer(const char opt[],const char text[],const char man[],PetscViewer *viewer,PetscViewerFormat *format,PetscBool  *set)
15538cc676e6SMatthew G Knepley {
15548cc676e6SMatthew G Knepley   PetscErrorCode ierr;
15558cc676e6SMatthew G Knepley   PetscOptions   amsopt;
15568cc676e6SMatthew G Knepley 
15578cc676e6SMatthew G Knepley   PetscFunctionBegin;
15588cc676e6SMatthew G Knepley   if (!PetscOptionsPublishCount) {
15598cc676e6SMatthew G Knepley     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
156064facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
15615b02f95dSBarry Smith     ierr = PetscStrdup("",(char**)&amsopt->data);CHKERRQ(ierr);
15628cc676e6SMatthew G Knepley   }
1563cffb1e40SBarry Smith   ierr = PetscOptionsGetViewer(PetscOptionsObject.comm,PetscOptionsObject.prefix,opt,viewer,format,set);CHKERRQ(ierr);
15648cc676e6SMatthew G Knepley   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
15658cc676e6SMatthew G Knepley     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%s>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,"",text,ManSection(man));CHKERRQ(ierr);
15668cc676e6SMatthew G Knepley   }
15678cc676e6SMatthew G Knepley   PetscFunctionReturn(0);
15688cc676e6SMatthew G Knepley }
15698cc676e6SMatthew G Knepley 
157053acd3b1SBarry Smith 
157153acd3b1SBarry Smith #undef __FUNCT__
157253acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsHead"
157353acd3b1SBarry Smith /*@C
1574b52f573bSBarry Smith      PetscOptionsHead - Puts a heading before listing any more published options. Used, for example,
157553acd3b1SBarry Smith             in KSPSetFromOptions_GMRES().
157653acd3b1SBarry Smith 
15773f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
157853acd3b1SBarry Smith 
157953acd3b1SBarry Smith    Input Parameter:
158053acd3b1SBarry Smith .   head - the heading text
158153acd3b1SBarry Smith 
158253acd3b1SBarry Smith 
158353acd3b1SBarry Smith    Level: intermediate
158453acd3b1SBarry Smith 
158553acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
158653acd3b1SBarry Smith 
1587b52f573bSBarry Smith           Can be followed by a call to PetscOptionsTail() in the same function.
158853acd3b1SBarry Smith 
158953acd3b1SBarry Smith    Concepts: options database^subheading
159053acd3b1SBarry Smith 
159153acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1592acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
159353acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
159453acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1595acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1596a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
159753acd3b1SBarry Smith @*/
15987087cfbeSBarry Smith PetscErrorCode  PetscOptionsHead(const char head[])
159953acd3b1SBarry Smith {
160053acd3b1SBarry Smith   PetscErrorCode ierr;
160153acd3b1SBarry Smith 
160253acd3b1SBarry Smith   PetscFunctionBegin;
160361b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
160453acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  %s\n",head);CHKERRQ(ierr);
160553acd3b1SBarry Smith   }
160653acd3b1SBarry Smith   PetscFunctionReturn(0);
160753acd3b1SBarry Smith }
160853acd3b1SBarry Smith 
160953acd3b1SBarry Smith 
161053acd3b1SBarry Smith 
161153acd3b1SBarry Smith 
161253acd3b1SBarry Smith 
161353acd3b1SBarry Smith 
1614