xref: /petsc/src/sys/objects/aoptions.c (revision 64bbc9efa26b1db25de9a163eb8454b1e25a5dcd)
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);
1635b02f95dSBarry Smith     tmp = (void*) 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 
376*64bbc9efSBarry Smith static const char *OptionsHeader = "<head>\n"
377*64bbc9efSBarry Smith                                    "<script jowererw type=\"text/javascript\" src=\"https://bitbucket.org/saws/saws/raw/master/js/jquery-1.9.1.js\"></script>\n"
378*64bbc9efSBarry Smith                                    "<script type=\"text/javascript\" src=\"https://bitbucket.org/saws/saws/raw/master/js/jsSAWs.js\"></script>\n"
379*64bbc9efSBarry Smith                                    "<script type=\"text/javascript\" src=\"PETSc.js\"></script>\n"
380*64bbc9efSBarry Smith                                    "<script>\n"
381*64bbc9efSBarry Smith                                       "jQuery(document).ready(function() {\n"
382*64bbc9efSBarry Smith                                          "PETSc.getAndDisplayDirectory(null,\"#variablesInfo\")\n"
383*64bbc9efSBarry Smith                                       "})\n"
384*64bbc9efSBarry Smith                                   "</script>\n"
385*64bbc9efSBarry Smith                                   "</head>\n";
386*64bbc9efSBarry Smith static const char *OptionsBodyBottom = "<div id=\"variablesInfo\" style=\"float:left\"></div>\n<br>\n</body>";
387*64bbc9efSBarry 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;
40788a9752cSBarry Smith   char           manname[16],textname[16];
4082657e9d9SBarry Smith   char           dir[1024];
409b3506946SBarry Smith 
410b3506946SBarry Smith   /* the next line is a bug, this will only work if all processors are here, the comm passed in is ignored!!! */
411b3506946SBarry Smith   sprintf(options,"Options_%d",count++);
412a297a907SKarl Rupp 
413e04113cfSBarry Smith   PetscOptionsObject.pprefix = PetscOptionsObject.prefix; /* SAWs will change this, so cannot pass prefix directly */
4141bc75a8dSBarry Smith 
415d50831c4SBarry Smith   ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s","_title");CHKERRQ(ierr);
4167781c08eSBarry Smith   PetscStackCallSAWs(SAWs_Register,(dir,&PetscOptionsObject.title,1,SAWs_READ,SAWs_STRING));
4177781c08eSBarry Smith   ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s","prefix");CHKERRQ(ierr);
4182657e9d9SBarry Smith   PetscStackCallSAWs(SAWs_Register,(dir,&PetscOptionsObject.pprefix,1,SAWs_READ,SAWs_STRING));
4192657e9d9SBarry Smith   PetscStackCallSAWs(SAWs_Register,("/PETSc/Options/ChangedMethod",&changedmethod,1,SAWs_WRITE,SAWs_BOOLEAN));
420b3506946SBarry Smith 
421b3506946SBarry Smith   while (next) {
422d50831c4SBarry Smith     sprintf(manname,"_man_%d",mancount);
4232657e9d9SBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",manname);CHKERRQ(ierr);
4247781c08eSBarry Smith     PetscStackCallSAWs(SAWs_Register,(dir,&next->man,1,SAWs_READ,SAWs_STRING));
425d50831c4SBarry Smith     sprintf(textname,"_text_%d",mancount++);
4267781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",textname);CHKERRQ(ierr);
4277781c08eSBarry Smith     PetscStackCallSAWs(SAWs_Register,(dir,&next->text,1,SAWs_READ,SAWs_STRING));
4289f32e415SBarry Smith 
429b3506946SBarry Smith     switch (next->type) {
430b3506946SBarry Smith     case OPTION_HEAD:
431b3506946SBarry Smith       break;
432b3506946SBarry Smith     case OPTION_INT_ARRAY:
4337781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4342657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_INT));
435b3506946SBarry Smith       break;
436b3506946SBarry Smith     case OPTION_REAL_ARRAY:
4377781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4382657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_DOUBLE));
439b3506946SBarry Smith       break;
440b3506946SBarry Smith     case OPTION_INT:
4417781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4422657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_INT));
443b3506946SBarry Smith       break;
444b3506946SBarry Smith     case OPTION_REAL:
4457781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4462657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_DOUBLE));
447b3506946SBarry Smith       break;
4487781c08eSBarry Smith     case OPTION_BOOL:
4497781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4502657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_BOOLEAN));
4511ae3d29cSBarry Smith       break;
4527781c08eSBarry Smith     case OPTION_BOOL_ARRAY:
4537781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4542657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_BOOLEAN));
45571f08665SBarry Smith       break;
456b3506946SBarry Smith     case OPTION_STRING:
4577781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4587781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
4591ae3d29cSBarry Smith       break;
4601ae3d29cSBarry Smith     case OPTION_STRING_ARRAY:
4617781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4622657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_STRING));
463b3506946SBarry Smith       break;
464a264d7a6SBarry Smith     case OPTION_FLIST:
465a264d7a6SBarry Smith       {
466a264d7a6SBarry Smith       PetscInt ntext;
4677781c08eSBarry Smith       ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4687781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
469a264d7a6SBarry Smith       ierr = PetscFunctionListGet(next->flist,(const char***)&next->edata,&ntext);CHKERRQ(ierr);
470a264d7a6SBarry Smith       PetscStackCallSAWs(SAWs_Set_Legal_Variable_Values,(dir,ntext,next->edata));
471a264d7a6SBarry Smith       }
472a264d7a6SBarry Smith       break;
4731ae3d29cSBarry Smith     case OPTION_ELIST:
474a264d7a6SBarry Smith       {
475a264d7a6SBarry Smith       PetscInt ntext = next->nlist;
4767781c08eSBarry Smith       ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4777781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
478ead66b60SBarry Smith       ierr = PetscMalloc1((ntext+1),(char***)&next->edata);CHKERRQ(ierr);
4791ae3d29cSBarry Smith       ierr = PetscMemcpy(next->edata,next->list,ntext*sizeof(char*));CHKERRQ(ierr);
480a264d7a6SBarry Smith       PetscStackCallSAWs(SAWs_Set_Legal_Variable_Values,(dir,ntext,next->edata));
481a264d7a6SBarry Smith       }
482a264d7a6SBarry Smith       break;
483b3506946SBarry Smith     default:
484b3506946SBarry Smith       break;
485b3506946SBarry Smith     }
486b3506946SBarry Smith     next = next->next;
487b3506946SBarry Smith   }
488b3506946SBarry Smith 
489b3506946SBarry Smith   /* wait until accessor has unlocked the memory */
490*64bbc9efSBarry Smith   PetscStackCallSAWs(SAWs_Push_Header,("index.html",OptionsHeader));
491*64bbc9efSBarry Smith   PetscStackCallSAWs(SAWs_Push_Body,("index.html",2,OptionsBodyBottom));
4927aab2a10SBarry Smith   ierr = PetscSAWsBlock();CHKERRQ(ierr);
493*64bbc9efSBarry Smith   PetscStackCallSAWs(SAWs_Pop_Header,("index.html"));
494*64bbc9efSBarry Smith   PetscStackCallSAWs(SAWs_Pop_Body,("index.html",2));
495b3506946SBarry Smith 
49688a9752cSBarry Smith   /* determine if any values have been set in GUI */
49788a9752cSBarry Smith   next = PetscOptionsObject.next;
49888a9752cSBarry Smith   while (next) {
49988a9752cSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
50088a9752cSBarry Smith     PetscStackCallSAWs(SAWs_Selected,(dir,&next->set));
50188a9752cSBarry Smith     next = next->next;
50288a9752cSBarry Smith   }
50388a9752cSBarry Smith 
504b3506946SBarry Smith   /* reset counter to -2; this updates the screen with the new options for the selected method */
505b3506946SBarry Smith   if (changedmethod) PetscOptionsPublishCount = -2;
506b3506946SBarry Smith 
5079a492a5cSBarry Smith   PetscStackCallSAWs(SAWs_Delete,("/PETSc/Options"));
508b3506946SBarry Smith   PetscFunctionReturn(0);
509b3506946SBarry Smith }
510b3506946SBarry Smith #endif
511b3506946SBarry Smith 
5126356e834SBarry Smith #undef __FUNCT__
51353acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEnd_Private"
51453acd3b1SBarry Smith PetscErrorCode PetscOptionsEnd_Private(void)
51553acd3b1SBarry Smith {
51653acd3b1SBarry Smith   PetscErrorCode ierr;
5176356e834SBarry Smith   PetscOptions   last;
5186356e834SBarry Smith   char           option[256],value[1024],tmp[32];
519330cf3c9SBarry Smith   size_t         j;
52053acd3b1SBarry Smith 
52153acd3b1SBarry Smith   PetscFunctionBegin;
522aee2cecaSBarry Smith   if (PetscOptionsObject.next) {
523b3506946SBarry Smith     if (!PetscOptionsPublishCount) {
524a264d7a6SBarry Smith #if defined(PETSC_HAVE_SAWS)
525475446a1SBarry Smith       ierr = PetscOptionsSAWsInput();CHKERRQ(ierr);
526b3506946SBarry Smith #else
52771f08665SBarry Smith       ierr = PetscOptionsGetFromTextInput();CHKERRQ(ierr);
528b3506946SBarry Smith #endif
529aee2cecaSBarry Smith     }
530aee2cecaSBarry Smith   }
5316356e834SBarry Smith 
532c31cb41cSBarry Smith   ierr = PetscFree(PetscOptionsObject.title);CHKERRQ(ierr);
533c31cb41cSBarry Smith   ierr = PetscFree(PetscOptionsObject.prefix);CHKERRQ(ierr);
5346356e834SBarry Smith 
5356356e834SBarry Smith   /* reset counter to -2; this updates the screen with the new options for the selected method */
5366356e834SBarry Smith   if (PetscOptionsObject.changedmethod) PetscOptionsPublishCount = -2;
53761b37b28SSatish Balay   /* reset alreadyprinted flag */
53861b37b28SSatish Balay   PetscOptionsObject.alreadyprinted = PETSC_FALSE;
5393194b578SJed Brown   if (PetscOptionsObject.object) PetscOptionsObject.object->optionsprinted = PETSC_TRUE;
5400298fd71SBarry Smith   PetscOptionsObject.object = NULL;
5416356e834SBarry Smith 
5426356e834SBarry Smith   while (PetscOptionsObject.next) {
5436356e834SBarry Smith     if (PetscOptionsObject.next->set) {
5446356e834SBarry Smith       if (PetscOptionsObject.prefix) {
5456356e834SBarry Smith         ierr = PetscStrcpy(option,"-");CHKERRQ(ierr);
5466356e834SBarry Smith         ierr = PetscStrcat(option,PetscOptionsObject.prefix);CHKERRQ(ierr);
5476356e834SBarry Smith         ierr = PetscStrcat(option,PetscOptionsObject.next->option+1);CHKERRQ(ierr);
5486356e834SBarry Smith       } else {
5496356e834SBarry Smith         ierr = PetscStrcpy(option,PetscOptionsObject.next->option);CHKERRQ(ierr);
5506356e834SBarry Smith       }
5516356e834SBarry Smith 
5526356e834SBarry Smith       switch (PetscOptionsObject.next->type) {
5536356e834SBarry Smith       case OPTION_HEAD:
5546356e834SBarry Smith         break;
555e26ddf31SBarry Smith       case OPTION_INT_ARRAY:
556e26ddf31SBarry Smith         sprintf(value,"%d",(int)((PetscInt*)PetscOptionsObject.next->data)[0]);
557e26ddf31SBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
558e26ddf31SBarry Smith           sprintf(tmp,"%d",(int)((PetscInt*)PetscOptionsObject.next->data)[j]);
559e26ddf31SBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
560e26ddf31SBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
561e26ddf31SBarry Smith         }
562e26ddf31SBarry Smith         break;
5636356e834SBarry Smith       case OPTION_INT:
5647a72a596SBarry Smith         sprintf(value,"%d",(int) *(PetscInt*)PetscOptionsObject.next->data);
5656356e834SBarry Smith         break;
5666356e834SBarry Smith       case OPTION_REAL:
5677a72a596SBarry Smith         sprintf(value,"%g",(double) *(PetscReal*)PetscOptionsObject.next->data);
5686356e834SBarry Smith         break;
5696356e834SBarry Smith       case OPTION_REAL_ARRAY:
5707a72a596SBarry Smith         sprintf(value,"%g",(double)((PetscReal*)PetscOptionsObject.next->data)[0]);
5716356e834SBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
5727a72a596SBarry Smith           sprintf(tmp,"%g",(double)((PetscReal*)PetscOptionsObject.next->data)[j]);
5736356e834SBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5746356e834SBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5756356e834SBarry Smith         }
5766356e834SBarry Smith         break;
5777781c08eSBarry Smith       case OPTION_BOOL:
57871f08665SBarry Smith         sprintf(value,"%d",*(int*)PetscOptionsObject.next->data);
5796356e834SBarry Smith         break;
5807781c08eSBarry Smith       case OPTION_BOOL_ARRAY:
581ace3abfcSBarry Smith         sprintf(value,"%d",(int)((PetscBool*)PetscOptionsObject.next->data)[0]);
5821ae3d29cSBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
583ace3abfcSBarry Smith           sprintf(tmp,"%d",(int)((PetscBool*)PetscOptionsObject.next->data)[j]);
5841ae3d29cSBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5851ae3d29cSBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5861ae3d29cSBarry Smith         }
5871ae3d29cSBarry Smith         break;
588a264d7a6SBarry Smith       case OPTION_FLIST:
5891ae3d29cSBarry Smith       case OPTION_ELIST:
5907781c08eSBarry Smith         ierr = PetscStrcpy(value,(char*)PetscOptionsObject.next->data);CHKERRQ(ierr);
5916356e834SBarry Smith         break;
5921ae3d29cSBarry Smith       case OPTION_STRING:
593475446a1SBarry Smith         ierr = PetscStrcpy(value,(char*)PetscOptionsObject.next->data);CHKERRQ(ierr);
5948da2146eSBarry Smith         break;
5951ae3d29cSBarry Smith       case OPTION_STRING_ARRAY:
5961ae3d29cSBarry Smith         sprintf(value,"%s",((char**)PetscOptionsObject.next->data)[0]);
5971ae3d29cSBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
5981ae3d29cSBarry Smith           sprintf(tmp,"%s",((char**)PetscOptionsObject.next->data)[j]);
5991ae3d29cSBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
6001ae3d29cSBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
6011ae3d29cSBarry Smith         }
6026356e834SBarry Smith         break;
6036356e834SBarry Smith       }
6046356e834SBarry Smith       ierr = PetscOptionsSetValue(option,value);CHKERRQ(ierr);
6056356e834SBarry Smith     }
606503cfb0cSBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->text);CHKERRQ(ierr);
607503cfb0cSBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->option);CHKERRQ(ierr);
6086356e834SBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->man);CHKERRQ(ierr);
60971f08665SBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->edata);CHKERRQ(ierr);
610c979a496SBarry Smith 
611c979a496SBarry Smith     if ((PetscOptionsObject.next->type == OPTION_STRING) || (PetscOptionsObject.next->type == OPTION_FLIST) || (PetscOptionsObject.next->type == OPTION_ELIST)){
61264facd6cSBarry Smith       /* must use system free since SAWs may have allocated it */
613c979a496SBarry Smith       free(PetscOptionsObject.next->data);
614c979a496SBarry Smith     } else {
6157781c08eSBarry Smith       ierr   = PetscFree(PetscOptionsObject.next->data);CHKERRQ(ierr);
616c979a496SBarry Smith     }
6177781c08eSBarry Smith 
6186356e834SBarry Smith     last                    = PetscOptionsObject.next;
6196356e834SBarry Smith     PetscOptionsObject.next = PetscOptionsObject.next->next;
6206356e834SBarry Smith     ierr                    = PetscFree(last);CHKERRQ(ierr);
6216356e834SBarry Smith   }
6226356e834SBarry Smith   PetscOptionsObject.next = 0;
62353acd3b1SBarry Smith   PetscFunctionReturn(0);
62453acd3b1SBarry Smith }
62553acd3b1SBarry Smith 
62653acd3b1SBarry Smith #undef __FUNCT__
62753acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEnum"
62853acd3b1SBarry Smith /*@C
62953acd3b1SBarry Smith    PetscOptionsEnum - Gets the enum value for a particular option in the database.
63053acd3b1SBarry Smith 
6313f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
63253acd3b1SBarry Smith 
63353acd3b1SBarry Smith    Input Parameters:
63453acd3b1SBarry Smith +  opt - option name
63553acd3b1SBarry Smith .  text - short string that describes the option
63653acd3b1SBarry Smith .  man - manual page with additional information on option
63753acd3b1SBarry Smith .  list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null
63853acd3b1SBarry Smith -  defaultv - the default (current) value
63953acd3b1SBarry Smith 
64053acd3b1SBarry Smith    Output Parameter:
64153acd3b1SBarry Smith +  value - the  value to return
642b32e0204SMatthew G Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
64353acd3b1SBarry Smith 
64453acd3b1SBarry Smith    Level: beginner
64553acd3b1SBarry Smith 
64653acd3b1SBarry Smith    Concepts: options database
64753acd3b1SBarry Smith 
64853acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
64953acd3b1SBarry Smith 
65053acd3b1SBarry Smith           list is usually something like PCASMTypes or some other predefined list of enum names
65153acd3b1SBarry Smith 
65253acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
653acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
654acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
65553acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
65653acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
657acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
658a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
65953acd3b1SBarry Smith @*/
6607087cfbeSBarry Smith PetscErrorCode  PetscOptionsEnum(const char opt[],const char text[],const char man[],const char *const *list,PetscEnum defaultv,PetscEnum *value,PetscBool  *set)
66153acd3b1SBarry Smith {
66253acd3b1SBarry Smith   PetscErrorCode ierr;
66353acd3b1SBarry Smith   PetscInt       ntext = 0;
664aa5bb8c0SSatish Balay   PetscInt       tval;
665ace3abfcSBarry Smith   PetscBool      tflg;
66653acd3b1SBarry Smith 
66753acd3b1SBarry Smith   PetscFunctionBegin;
66853acd3b1SBarry Smith   while (list[ntext++]) {
669e32f2f54SBarry Smith     if (ntext > 50) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument appears to be wrong or have more than 50 entries");
67053acd3b1SBarry Smith   }
671e32f2f54SBarry Smith   if (ntext < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument must have at least two entries: typename and type prefix");
67253acd3b1SBarry Smith   ntext -= 3;
673aa5bb8c0SSatish Balay   ierr   = PetscOptionsEList(opt,text,man,list,ntext,list[defaultv],&tval,&tflg);CHKERRQ(ierr);
674aa5bb8c0SSatish Balay   /* with PETSC_USE_64BIT_INDICES sizeof(PetscInt) != sizeof(PetscEnum) */
675aa5bb8c0SSatish Balay   if (tflg) *value = (PetscEnum)tval;
676aa5bb8c0SSatish Balay   if (set)  *set   = tflg;
67753acd3b1SBarry Smith   PetscFunctionReturn(0);
67853acd3b1SBarry Smith }
67953acd3b1SBarry Smith 
68053acd3b1SBarry Smith /* -------------------------------------------------------------------------------------------------------------*/
68153acd3b1SBarry Smith #undef __FUNCT__
68253acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsInt"
68353acd3b1SBarry Smith /*@C
68453acd3b1SBarry Smith    PetscOptionsInt - Gets the integer value for a particular option in the database.
68553acd3b1SBarry Smith 
6863f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
68753acd3b1SBarry Smith 
68853acd3b1SBarry Smith    Input Parameters:
68953acd3b1SBarry Smith +  opt - option name
69053acd3b1SBarry Smith .  text - short string that describes the option
69153acd3b1SBarry Smith .  man - manual page with additional information on option
69253acd3b1SBarry Smith -  defaultv - the default (current) value
69353acd3b1SBarry Smith 
69453acd3b1SBarry Smith    Output Parameter:
69553acd3b1SBarry Smith +  value - the integer value to return
69653acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
69753acd3b1SBarry Smith 
69853acd3b1SBarry Smith    Level: beginner
69953acd3b1SBarry Smith 
70053acd3b1SBarry Smith    Concepts: options database^has int
70153acd3b1SBarry Smith 
70253acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
70353acd3b1SBarry Smith 
70453acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
705acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
706acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
70753acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
70853acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
709acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
710a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
71153acd3b1SBarry Smith @*/
7127087cfbeSBarry Smith PetscErrorCode  PetscOptionsInt(const char opt[],const char text[],const char man[],PetscInt defaultv,PetscInt *value,PetscBool  *set)
71353acd3b1SBarry Smith {
71453acd3b1SBarry Smith   PetscErrorCode ierr;
7156356e834SBarry Smith   PetscOptions   amsopt;
71653acd3b1SBarry Smith 
71753acd3b1SBarry Smith   PetscFunctionBegin;
718b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
7196356e834SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_INT,&amsopt);CHKERRQ(ierr);
7206356e834SBarry Smith     ierr = PetscMalloc(sizeof(PetscInt),&amsopt->data);CHKERRQ(ierr);
721a297a907SKarl Rupp 
7226356e834SBarry Smith     *(PetscInt*)amsopt->data = defaultv;
723af6d86caSBarry Smith   }
72453acd3b1SBarry Smith   ierr = PetscOptionsGetInt(PetscOptionsObject.prefix,opt,value,set);CHKERRQ(ierr);
72561b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
7262aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%d>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,defaultv,text,ManSection(man));CHKERRQ(ierr);
72753acd3b1SBarry Smith   }
72853acd3b1SBarry Smith   PetscFunctionReturn(0);
72953acd3b1SBarry Smith }
73053acd3b1SBarry Smith 
73153acd3b1SBarry Smith #undef __FUNCT__
73253acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsString"
73353acd3b1SBarry Smith /*@C
73453acd3b1SBarry Smith    PetscOptionsString - Gets the string value for a particular option in the database.
73553acd3b1SBarry Smith 
7363f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
73753acd3b1SBarry Smith 
73853acd3b1SBarry Smith    Input Parameters:
73953acd3b1SBarry Smith +  opt - option name
74053acd3b1SBarry Smith .  text - short string that describes the option
74153acd3b1SBarry Smith .  man - manual page with additional information on option
742bcbf2dc5SJed Brown .  defaultv - the default (current) value
743bcbf2dc5SJed Brown -  len - length of the result string including null terminator
74453acd3b1SBarry Smith 
74553acd3b1SBarry Smith    Output Parameter:
74653acd3b1SBarry Smith +  value - the value to return
74753acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
74853acd3b1SBarry Smith 
74953acd3b1SBarry Smith    Level: beginner
75053acd3b1SBarry Smith 
75153acd3b1SBarry Smith    Concepts: options database^has int
75253acd3b1SBarry Smith 
75353acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
75453acd3b1SBarry Smith 
7557fccdfe4SBarry 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).
7567fccdfe4SBarry Smith 
75753acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
758acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
759acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsReal(), PetscOptionsBool(),
76053acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
76153acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
762acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
763a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
76453acd3b1SBarry Smith @*/
7657087cfbeSBarry Smith PetscErrorCode  PetscOptionsString(const char opt[],const char text[],const char man[],const char defaultv[],char value[],size_t len,PetscBool  *set)
76653acd3b1SBarry Smith {
76753acd3b1SBarry Smith   PetscErrorCode ierr;
768aee2cecaSBarry Smith   PetscOptions   amsopt;
76953acd3b1SBarry Smith 
77053acd3b1SBarry Smith   PetscFunctionBegin;
771b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
772aee2cecaSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
77364facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
7745b02f95dSBarry Smith     ierr = PetscStrdup(defaultv ? defaultv : "",(char**)&amsopt->data);CHKERRQ(ierr);
775af6d86caSBarry Smith   }
77653acd3b1SBarry Smith   ierr = PetscOptionsGetString(PetscOptionsObject.prefix,opt,value,len,set);CHKERRQ(ierr);
77761b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
7782aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%s>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,defaultv,text,ManSection(man));CHKERRQ(ierr);
77953acd3b1SBarry Smith   }
78053acd3b1SBarry Smith   PetscFunctionReturn(0);
78153acd3b1SBarry Smith }
78253acd3b1SBarry Smith 
78353acd3b1SBarry Smith #undef __FUNCT__
78453acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsReal"
78553acd3b1SBarry Smith /*@C
78653acd3b1SBarry Smith    PetscOptionsReal - Gets the PetscReal value for a particular option in the database.
78753acd3b1SBarry Smith 
7883f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
78953acd3b1SBarry Smith 
79053acd3b1SBarry Smith    Input Parameters:
79153acd3b1SBarry Smith +  opt - option name
79253acd3b1SBarry Smith .  text - short string that describes the option
79353acd3b1SBarry Smith .  man - manual page with additional information on option
79453acd3b1SBarry Smith -  defaultv - the default (current) value
79553acd3b1SBarry Smith 
79653acd3b1SBarry Smith    Output Parameter:
79753acd3b1SBarry Smith +  value - the value to return
79853acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
79953acd3b1SBarry Smith 
80053acd3b1SBarry Smith    Level: beginner
80153acd3b1SBarry Smith 
80253acd3b1SBarry Smith    Concepts: options database^has int
80353acd3b1SBarry Smith 
80453acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
80553acd3b1SBarry Smith 
80653acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
807acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
808acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
80953acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
81053acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
811acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
812a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
81353acd3b1SBarry Smith @*/
8147087cfbeSBarry Smith PetscErrorCode  PetscOptionsReal(const char opt[],const char text[],const char man[],PetscReal defaultv,PetscReal *value,PetscBool  *set)
81553acd3b1SBarry Smith {
81653acd3b1SBarry Smith   PetscErrorCode ierr;
817538aa990SBarry Smith   PetscOptions   amsopt;
81853acd3b1SBarry Smith 
81953acd3b1SBarry Smith   PetscFunctionBegin;
820b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
821538aa990SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_REAL,&amsopt);CHKERRQ(ierr);
822538aa990SBarry Smith     ierr = PetscMalloc(sizeof(PetscReal),&amsopt->data);CHKERRQ(ierr);
823a297a907SKarl Rupp 
824538aa990SBarry Smith     *(PetscReal*)amsopt->data = defaultv;
825538aa990SBarry Smith   }
82653acd3b1SBarry Smith   ierr = PetscOptionsGetReal(PetscOptionsObject.prefix,opt,value,set);CHKERRQ(ierr);
82761b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
8282aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%G>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,defaultv,text,ManSection(man));CHKERRQ(ierr);
82953acd3b1SBarry Smith   }
83053acd3b1SBarry Smith   PetscFunctionReturn(0);
83153acd3b1SBarry Smith }
83253acd3b1SBarry Smith 
83353acd3b1SBarry Smith #undef __FUNCT__
83453acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsScalar"
83553acd3b1SBarry Smith /*@C
83653acd3b1SBarry Smith    PetscOptionsScalar - Gets the scalar value for a particular option in the database.
83753acd3b1SBarry Smith 
8383f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
83953acd3b1SBarry Smith 
84053acd3b1SBarry Smith    Input Parameters:
84153acd3b1SBarry Smith +  opt - option name
84253acd3b1SBarry Smith .  text - short string that describes the option
84353acd3b1SBarry Smith .  man - manual page with additional information on option
84453acd3b1SBarry Smith -  defaultv - the default (current) value
84553acd3b1SBarry Smith 
84653acd3b1SBarry Smith    Output Parameter:
84753acd3b1SBarry Smith +  value - the value to return
84853acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
84953acd3b1SBarry Smith 
85053acd3b1SBarry Smith    Level: beginner
85153acd3b1SBarry Smith 
85253acd3b1SBarry Smith    Concepts: options database^has int
85353acd3b1SBarry Smith 
85453acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
85553acd3b1SBarry Smith 
85653acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
857acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
858acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
85953acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
86053acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
861acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
862a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
86353acd3b1SBarry Smith @*/
8647087cfbeSBarry Smith PetscErrorCode  PetscOptionsScalar(const char opt[],const char text[],const char man[],PetscScalar defaultv,PetscScalar *value,PetscBool  *set)
86553acd3b1SBarry Smith {
86653acd3b1SBarry Smith   PetscErrorCode ierr;
86753acd3b1SBarry Smith 
86853acd3b1SBarry Smith   PetscFunctionBegin;
86953acd3b1SBarry Smith #if !defined(PETSC_USE_COMPLEX)
87053acd3b1SBarry Smith   ierr = PetscOptionsReal(opt,text,man,defaultv,value,set);CHKERRQ(ierr);
87153acd3b1SBarry Smith #else
87253acd3b1SBarry Smith   ierr = PetscOptionsGetScalar(PetscOptionsObject.prefix,opt,value,set);CHKERRQ(ierr);
87353acd3b1SBarry Smith #endif
87453acd3b1SBarry Smith   PetscFunctionReturn(0);
87553acd3b1SBarry Smith }
87653acd3b1SBarry Smith 
87753acd3b1SBarry Smith #undef __FUNCT__
87853acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsName"
87953acd3b1SBarry Smith /*@C
88090d69ab7SBarry 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
88190d69ab7SBarry Smith                       its value is set to false.
88253acd3b1SBarry Smith 
8833f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
88453acd3b1SBarry Smith 
88553acd3b1SBarry Smith    Input Parameters:
88653acd3b1SBarry Smith +  opt - option name
88753acd3b1SBarry Smith .  text - short string that describes the option
88853acd3b1SBarry Smith -  man - manual page with additional information on option
88953acd3b1SBarry Smith 
89053acd3b1SBarry Smith    Output Parameter:
89153acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
89253acd3b1SBarry Smith 
89353acd3b1SBarry Smith    Level: beginner
89453acd3b1SBarry Smith 
89553acd3b1SBarry Smith    Concepts: options database^has int
89653acd3b1SBarry Smith 
89753acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
89853acd3b1SBarry Smith 
89953acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
900acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
901acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
90253acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
90353acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
904acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
905a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
90653acd3b1SBarry Smith @*/
9077087cfbeSBarry Smith PetscErrorCode  PetscOptionsName(const char opt[],const char text[],const char man[],PetscBool  *flg)
90853acd3b1SBarry Smith {
90953acd3b1SBarry Smith   PetscErrorCode ierr;
9101ae3d29cSBarry Smith   PetscOptions   amsopt;
91153acd3b1SBarry Smith 
91253acd3b1SBarry Smith   PetscFunctionBegin;
9131ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
9147781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
915ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
916a297a907SKarl Rupp 
917ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
9181ae3d29cSBarry Smith   }
91953acd3b1SBarry Smith   ierr = PetscOptionsHasName(PetscOptionsObject.prefix,opt,flg);CHKERRQ(ierr);
92061b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
9212aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
92253acd3b1SBarry Smith   }
92353acd3b1SBarry Smith   PetscFunctionReturn(0);
92453acd3b1SBarry Smith }
92553acd3b1SBarry Smith 
92653acd3b1SBarry Smith #undef __FUNCT__
927a264d7a6SBarry Smith #define __FUNCT__ "PetscOptionsFList"
92853acd3b1SBarry Smith /*@C
929a264d7a6SBarry Smith      PetscOptionsFList - Puts a list of option values that a single one may be selected from
93053acd3b1SBarry Smith 
9313f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
93253acd3b1SBarry Smith 
93353acd3b1SBarry Smith    Input Parameters:
93453acd3b1SBarry Smith +  opt - option name
93553acd3b1SBarry Smith .  text - short string that describes the option
93653acd3b1SBarry Smith .  man - manual page with additional information on option
93753acd3b1SBarry Smith .  list - the possible choices
9383cc1e11dSBarry Smith .  defaultv - the default (current) value
9393cc1e11dSBarry Smith -  len - the length of the character array value
94053acd3b1SBarry Smith 
94153acd3b1SBarry Smith    Output Parameter:
94253acd3b1SBarry Smith +  value - the value to return
94353acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
94453acd3b1SBarry Smith 
94553acd3b1SBarry Smith    Level: intermediate
94653acd3b1SBarry Smith 
94753acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
94853acd3b1SBarry Smith 
94953acd3b1SBarry Smith    See PetscOptionsEList() for when the choices are given in a string array
95053acd3b1SBarry Smith 
95153acd3b1SBarry Smith    To get a listing of all currently specified options,
95288c29154SBarry Smith     see PetscOptionsView() or PetscOptionsGetAll()
95353acd3b1SBarry Smith 
954eabe10d7SBarry Smith    Developer Note: This cannot check for invalid selection because of things like MATAIJ that are not included in the list
955eabe10d7SBarry Smith 
95653acd3b1SBarry Smith    Concepts: options database^list
95753acd3b1SBarry Smith 
95853acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
959acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
96053acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
96153acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
962acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
963a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsEnum()
96453acd3b1SBarry Smith @*/
965a264d7a6SBarry Smith PetscErrorCode  PetscOptionsFList(const char opt[],const char ltext[],const char man[],PetscFunctionList list,const char defaultv[],char value[],size_t len,PetscBool  *set)
96653acd3b1SBarry Smith {
96753acd3b1SBarry Smith   PetscErrorCode ierr;
9683cc1e11dSBarry Smith   PetscOptions   amsopt;
96953acd3b1SBarry Smith 
97053acd3b1SBarry Smith   PetscFunctionBegin;
971b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
972a264d7a6SBarry Smith     ierr = PetscOptionsCreate_Private(opt,ltext,man,OPTION_FLIST,&amsopt);CHKERRQ(ierr);
97364facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
9745b02f95dSBarry Smith     ierr = PetscStrdup(defaultv ? defaultv : "",(char**)&amsopt->data);CHKERRQ(ierr);
9753cc1e11dSBarry Smith     amsopt->flist = list;
9763cc1e11dSBarry Smith   }
97753acd3b1SBarry Smith   ierr = PetscOptionsGetString(PetscOptionsObject.prefix,opt,value,len,set);CHKERRQ(ierr);
97861b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
979140e18c1SBarry Smith     ierr = PetscFunctionListPrintTypes(PetscOptionsObject.comm,stdout,PetscOptionsObject.prefix,opt,ltext,man,list,defaultv);CHKERRQ(ierr);CHKERRQ(ierr);
98053acd3b1SBarry Smith   }
98153acd3b1SBarry Smith   PetscFunctionReturn(0);
98253acd3b1SBarry Smith }
98353acd3b1SBarry Smith 
98453acd3b1SBarry Smith #undef __FUNCT__
98553acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEList"
98653acd3b1SBarry Smith /*@C
98753acd3b1SBarry Smith      PetscOptionsEList - Puts a list of option values that a single one may be selected from
98853acd3b1SBarry Smith 
9893f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
99053acd3b1SBarry Smith 
99153acd3b1SBarry Smith    Input Parameters:
99253acd3b1SBarry Smith +  opt - option name
99353acd3b1SBarry Smith .  ltext - short string that describes the option
99453acd3b1SBarry Smith .  man - manual page with additional information on option
995a264d7a6SBarry Smith .  list - the possible choices (one of these must be selected, anything else is invalid)
99653acd3b1SBarry Smith .  ntext - number of choices
99753acd3b1SBarry Smith -  defaultv - the default (current) value
99853acd3b1SBarry Smith 
99953acd3b1SBarry Smith    Output Parameter:
100053acd3b1SBarry Smith +  value - the index of the value to return
100153acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
100253acd3b1SBarry Smith 
100353acd3b1SBarry Smith    Level: intermediate
100453acd3b1SBarry Smith 
100553acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
100653acd3b1SBarry Smith 
1007a264d7a6SBarry Smith    See PetscOptionsFList() for when the choices are given in a PetscFunctionList()
100853acd3b1SBarry Smith 
100953acd3b1SBarry Smith    Concepts: options database^list
101053acd3b1SBarry Smith 
101153acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1012acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
101353acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
101453acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1015acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1016a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEnum()
101753acd3b1SBarry Smith @*/
10187087cfbeSBarry 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)
101953acd3b1SBarry Smith {
102053acd3b1SBarry Smith   PetscErrorCode ierr;
102153acd3b1SBarry Smith   PetscInt       i;
10221ae3d29cSBarry Smith   PetscOptions   amsopt;
102353acd3b1SBarry Smith 
102453acd3b1SBarry Smith   PetscFunctionBegin;
10251ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
1026d5649816SBarry Smith     ierr = PetscOptionsCreate_Private(opt,ltext,man,OPTION_ELIST,&amsopt);CHKERRQ(ierr);
102764facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
10285b02f95dSBarry Smith     ierr = PetscStrdup(defaultv ? defaultv : "",(char**)&amsopt->data);CHKERRQ(ierr);
10291ae3d29cSBarry Smith     amsopt->list  = list;
10301ae3d29cSBarry Smith     amsopt->nlist = ntext;
10311ae3d29cSBarry Smith   }
103253acd3b1SBarry Smith   ierr = PetscOptionsGetEList(PetscOptionsObject.prefix,opt,list,ntext,value,set);CHKERRQ(ierr);
103361b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
103453acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%s> (choose one of)",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,defaultv);CHKERRQ(ierr);
103553acd3b1SBarry Smith     for (i=0; i<ntext; i++) {
103653acd3b1SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm," %s",list[i]);CHKERRQ(ierr);
103753acd3b1SBarry Smith     }
10382aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm," (%s)\n",ManSection(man));CHKERRQ(ierr);
103953acd3b1SBarry Smith   }
104053acd3b1SBarry Smith   PetscFunctionReturn(0);
104153acd3b1SBarry Smith }
104253acd3b1SBarry Smith 
104353acd3b1SBarry Smith #undef __FUNCT__
1044acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolGroupBegin"
104553acd3b1SBarry Smith /*@C
1046acfcf0e5SJed Brown      PetscOptionsBoolGroupBegin - First in a series of logical queries on the options database for
1047d5649816SBarry Smith        which at most a single value can be true.
104853acd3b1SBarry Smith 
10493f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
105053acd3b1SBarry Smith 
105153acd3b1SBarry Smith    Input Parameters:
105253acd3b1SBarry Smith +  opt - option name
105353acd3b1SBarry Smith .  text - short string that describes the option
105453acd3b1SBarry Smith -  man - manual page with additional information on option
105553acd3b1SBarry Smith 
105653acd3b1SBarry Smith    Output Parameter:
105753acd3b1SBarry Smith .  flg - whether that option was set or not
105853acd3b1SBarry Smith 
105953acd3b1SBarry Smith    Level: intermediate
106053acd3b1SBarry Smith 
106153acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
106253acd3b1SBarry Smith 
1063acfcf0e5SJed Brown    Must be followed by 0 or more PetscOptionsBoolGroup()s and PetscOptionsBoolGroupEnd()
106453acd3b1SBarry Smith 
106553acd3b1SBarry Smith     Concepts: options database^logical group
106653acd3b1SBarry Smith 
106753acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1068acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
106953acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
107053acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1071acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1072a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
107353acd3b1SBarry Smith @*/
10747087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolGroupBegin(const char opt[],const char text[],const char man[],PetscBool  *flg)
107553acd3b1SBarry Smith {
107653acd3b1SBarry Smith   PetscErrorCode ierr;
10771ae3d29cSBarry Smith   PetscOptions   amsopt;
107853acd3b1SBarry Smith 
107953acd3b1SBarry Smith   PetscFunctionBegin;
10801ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
10817781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1082ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1083a297a907SKarl Rupp 
1084ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
10851ae3d29cSBarry Smith   }
108668b16fdaSBarry Smith   *flg = PETSC_FALSE;
10870298fd71SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,NULL);CHKERRQ(ierr);
108861b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
108953acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  Pick at most one of -------------\n");CHKERRQ(ierr);
10902aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"    -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
109153acd3b1SBarry Smith   }
109253acd3b1SBarry Smith   PetscFunctionReturn(0);
109353acd3b1SBarry Smith }
109453acd3b1SBarry Smith 
109553acd3b1SBarry Smith #undef __FUNCT__
1096acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolGroup"
109753acd3b1SBarry Smith /*@C
1098acfcf0e5SJed Brown      PetscOptionsBoolGroup - One in a series of logical queries on the options database for
1099d5649816SBarry Smith        which at most a single value can be true.
110053acd3b1SBarry Smith 
11013f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
110253acd3b1SBarry Smith 
110353acd3b1SBarry Smith    Input Parameters:
110453acd3b1SBarry Smith +  opt - option name
110553acd3b1SBarry Smith .  text - short string that describes the option
110653acd3b1SBarry Smith -  man - manual page with additional information on option
110753acd3b1SBarry Smith 
110853acd3b1SBarry Smith    Output Parameter:
110953acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
111053acd3b1SBarry Smith 
111153acd3b1SBarry Smith    Level: intermediate
111253acd3b1SBarry Smith 
111353acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
111453acd3b1SBarry Smith 
1115acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin() and preceded a PetscOptionsBoolGroupEnd()
111653acd3b1SBarry Smith 
111753acd3b1SBarry Smith     Concepts: options database^logical group
111853acd3b1SBarry Smith 
111953acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1120acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
112153acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
112253acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1123acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1124a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
112553acd3b1SBarry Smith @*/
11267087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolGroup(const char opt[],const char text[],const char man[],PetscBool  *flg)
112753acd3b1SBarry Smith {
112853acd3b1SBarry Smith   PetscErrorCode ierr;
11291ae3d29cSBarry Smith   PetscOptions   amsopt;
113053acd3b1SBarry Smith 
113153acd3b1SBarry Smith   PetscFunctionBegin;
11321ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
11337781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1134ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1135a297a907SKarl Rupp 
1136ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
11371ae3d29cSBarry Smith   }
113817326d04SJed Brown   *flg = PETSC_FALSE;
11390298fd71SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,NULL);CHKERRQ(ierr);
114061b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
11412aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"    -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
114253acd3b1SBarry Smith   }
114353acd3b1SBarry Smith   PetscFunctionReturn(0);
114453acd3b1SBarry Smith }
114553acd3b1SBarry Smith 
114653acd3b1SBarry Smith #undef __FUNCT__
1147acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolGroupEnd"
114853acd3b1SBarry Smith /*@C
1149acfcf0e5SJed Brown      PetscOptionsBoolGroupEnd - Last in a series of logical queries on the options database for
1150d5649816SBarry Smith        which at most a single value can be true.
115153acd3b1SBarry Smith 
11523f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
115353acd3b1SBarry Smith 
115453acd3b1SBarry Smith    Input Parameters:
115553acd3b1SBarry Smith +  opt - option name
115653acd3b1SBarry Smith .  text - short string that describes the option
115753acd3b1SBarry Smith -  man - manual page with additional information on option
115853acd3b1SBarry Smith 
115953acd3b1SBarry Smith    Output Parameter:
116053acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
116153acd3b1SBarry Smith 
116253acd3b1SBarry Smith    Level: intermediate
116353acd3b1SBarry Smith 
116453acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
116553acd3b1SBarry Smith 
1166acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin()
116753acd3b1SBarry Smith 
116853acd3b1SBarry Smith     Concepts: options database^logical group
116953acd3b1SBarry Smith 
117053acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1171acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
117253acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
117353acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1174acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1175a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
117653acd3b1SBarry Smith @*/
11777087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolGroupEnd(const char opt[],const char text[],const char man[],PetscBool  *flg)
117853acd3b1SBarry Smith {
117953acd3b1SBarry Smith   PetscErrorCode ierr;
11801ae3d29cSBarry Smith   PetscOptions   amsopt;
118153acd3b1SBarry Smith 
118253acd3b1SBarry Smith   PetscFunctionBegin;
11831ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
11847781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1185ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1186a297a907SKarl Rupp 
1187ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
11881ae3d29cSBarry Smith   }
118917326d04SJed Brown   *flg = PETSC_FALSE;
11900298fd71SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,NULL);CHKERRQ(ierr);
119161b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
11922aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"    -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
119353acd3b1SBarry Smith   }
119453acd3b1SBarry Smith   PetscFunctionReturn(0);
119553acd3b1SBarry Smith }
119653acd3b1SBarry Smith 
119753acd3b1SBarry Smith #undef __FUNCT__
1198acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBool"
119953acd3b1SBarry Smith /*@C
1200acfcf0e5SJed Brown    PetscOptionsBool - Determines if a particular option is in the database with a true or false
120153acd3b1SBarry Smith 
12023f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
120353acd3b1SBarry Smith 
120453acd3b1SBarry Smith    Input Parameters:
120553acd3b1SBarry Smith +  opt - option name
120653acd3b1SBarry Smith .  text - short string that describes the option
120753acd3b1SBarry Smith -  man - manual page with additional information on option
120853acd3b1SBarry Smith 
120953acd3b1SBarry Smith    Output Parameter:
121053acd3b1SBarry Smith .  flg - PETSC_TRUE or PETSC_FALSE
121153acd3b1SBarry Smith .  set - PETSC_TRUE if found, else PETSC_FALSE
121253acd3b1SBarry Smith 
121353acd3b1SBarry Smith    Level: beginner
121453acd3b1SBarry Smith 
121553acd3b1SBarry Smith    Concepts: options database^logical
121653acd3b1SBarry Smith 
121753acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
121853acd3b1SBarry Smith 
121953acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
1220acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
1221acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
122253acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
122353acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1224acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1225a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
122653acd3b1SBarry Smith @*/
12277087cfbeSBarry Smith PetscErrorCode  PetscOptionsBool(const char opt[],const char text[],const char man[],PetscBool deflt,PetscBool  *flg,PetscBool  *set)
122853acd3b1SBarry Smith {
122953acd3b1SBarry Smith   PetscErrorCode ierr;
1230ace3abfcSBarry Smith   PetscBool      iset;
1231aee2cecaSBarry Smith   PetscOptions   amsopt;
123253acd3b1SBarry Smith 
123353acd3b1SBarry Smith   PetscFunctionBegin;
1234b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
12357781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1236ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1237a297a907SKarl Rupp 
1238ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = deflt;
1239af6d86caSBarry Smith   }
1240acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,&iset);CHKERRQ(ierr);
124153acd3b1SBarry Smith   if (!iset) {
124253acd3b1SBarry Smith     if (flg) *flg = deflt;
124353acd3b1SBarry Smith   }
124453acd3b1SBarry Smith   if (set) *set = iset;
124561b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
1246ace3abfcSBarry Smith     const char *v = PetscBools[deflt];
12472aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s: <%s> %s (%s)\n",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,v,text,ManSection(man));CHKERRQ(ierr);
124853acd3b1SBarry Smith   }
124953acd3b1SBarry Smith   PetscFunctionReturn(0);
125053acd3b1SBarry Smith }
125153acd3b1SBarry Smith 
125253acd3b1SBarry Smith #undef __FUNCT__
125353acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsRealArray"
125453acd3b1SBarry Smith /*@C
125553acd3b1SBarry Smith    PetscOptionsRealArray - Gets an array of double values for a particular
125653acd3b1SBarry Smith    option in the database. The values must be separated with commas with
125753acd3b1SBarry Smith    no intervening spaces.
125853acd3b1SBarry Smith 
12593f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
126053acd3b1SBarry Smith 
126153acd3b1SBarry Smith    Input Parameters:
126253acd3b1SBarry Smith +  opt - the option one is seeking
126353acd3b1SBarry Smith .  text - short string describing option
126453acd3b1SBarry Smith .  man - manual page for option
126553acd3b1SBarry Smith -  nmax - maximum number of values
126653acd3b1SBarry Smith 
126753acd3b1SBarry Smith    Output Parameter:
126853acd3b1SBarry Smith +  value - location to copy values
126953acd3b1SBarry Smith .  nmax - actual number of values found
127053acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
127153acd3b1SBarry Smith 
127253acd3b1SBarry Smith    Level: beginner
127353acd3b1SBarry Smith 
127453acd3b1SBarry Smith    Notes:
127553acd3b1SBarry Smith    The user should pass in an array of doubles
127653acd3b1SBarry Smith 
127753acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
127853acd3b1SBarry Smith 
127953acd3b1SBarry Smith    Concepts: options database^array of strings
128053acd3b1SBarry Smith 
128153acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1282acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
128353acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
128453acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1285acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1286a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
128753acd3b1SBarry Smith @*/
12887087cfbeSBarry Smith PetscErrorCode  PetscOptionsRealArray(const char opt[],const char text[],const char man[],PetscReal value[],PetscInt *n,PetscBool  *set)
128953acd3b1SBarry Smith {
129053acd3b1SBarry Smith   PetscErrorCode ierr;
129153acd3b1SBarry Smith   PetscInt       i;
1292e26ddf31SBarry Smith   PetscOptions   amsopt;
129353acd3b1SBarry Smith 
129453acd3b1SBarry Smith   PetscFunctionBegin;
1295b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
1296e26ddf31SBarry Smith     PetscReal *vals;
1297e26ddf31SBarry Smith 
1298e26ddf31SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_REAL_ARRAY,&amsopt);CHKERRQ(ierr);
12996e02237eSPeter Brune     ierr = PetscMalloc1((*n),(PetscReal**)&amsopt->data);CHKERRQ(ierr);
1300e26ddf31SBarry Smith     vals = (PetscReal*)amsopt->data;
1301e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1302e26ddf31SBarry Smith     amsopt->arraylength = *n;
1303e26ddf31SBarry Smith   }
130453acd3b1SBarry Smith   ierr = PetscOptionsGetRealArray(PetscOptionsObject.prefix,opt,value,n,set);CHKERRQ(ierr);
130561b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
1306a83599f4SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%G",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,value[0]);CHKERRQ(ierr);
130753acd3b1SBarry Smith     for (i=1; i<*n; i++) {
1308a83599f4SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,",%G",value[i]);CHKERRQ(ierr);
130953acd3b1SBarry Smith     }
13102aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
131153acd3b1SBarry Smith   }
131253acd3b1SBarry Smith   PetscFunctionReturn(0);
131353acd3b1SBarry Smith }
131453acd3b1SBarry Smith 
131553acd3b1SBarry Smith 
131653acd3b1SBarry Smith #undef __FUNCT__
131753acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsIntArray"
131853acd3b1SBarry Smith /*@C
131953acd3b1SBarry Smith    PetscOptionsIntArray - Gets an array of integers for a particular
1320b32a342fSShri Abhyankar    option in the database.
132153acd3b1SBarry Smith 
13223f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
132353acd3b1SBarry Smith 
132453acd3b1SBarry Smith    Input Parameters:
132553acd3b1SBarry Smith +  opt - the option one is seeking
132653acd3b1SBarry Smith .  text - short string describing option
132753acd3b1SBarry Smith .  man - manual page for option
1328f8a50e2bSBarry Smith -  n - maximum number of values
132953acd3b1SBarry Smith 
133053acd3b1SBarry Smith    Output Parameter:
133153acd3b1SBarry Smith +  value - location to copy values
1332f8a50e2bSBarry Smith .  n - actual number of values found
133353acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
133453acd3b1SBarry Smith 
133553acd3b1SBarry Smith    Level: beginner
133653acd3b1SBarry Smith 
133753acd3b1SBarry Smith    Notes:
1338b32a342fSShri Abhyankar    The array can be passed as
1339b32a342fSShri Abhyankar    a comma seperated list:                                 0,1,2,3,4,5,6,7
13400fd488f5SShri Abhyankar    a range (start-end+1):                                  0-8
13410fd488f5SShri Abhyankar    a range with given increment (start-end+1:inc):         0-7:2
13420fd488f5SShri Abhyankar    a combination of values and ranges seperated by commas: 0,1-8,8-15:2
1343b32a342fSShri Abhyankar 
1344b32a342fSShri Abhyankar    There must be no intervening spaces between the values.
134553acd3b1SBarry Smith 
134653acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
134753acd3b1SBarry Smith 
1348b32a342fSShri Abhyankar    Concepts: options database^array of ints
134953acd3b1SBarry Smith 
135053acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1351acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
135253acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
135353acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1354acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1355a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsRealArray()
135653acd3b1SBarry Smith @*/
13577087cfbeSBarry Smith PetscErrorCode  PetscOptionsIntArray(const char opt[],const char text[],const char man[],PetscInt value[],PetscInt *n,PetscBool  *set)
135853acd3b1SBarry Smith {
135953acd3b1SBarry Smith   PetscErrorCode ierr;
136053acd3b1SBarry Smith   PetscInt       i;
1361e26ddf31SBarry Smith   PetscOptions   amsopt;
136253acd3b1SBarry Smith 
136353acd3b1SBarry Smith   PetscFunctionBegin;
1364b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
1365e26ddf31SBarry Smith     PetscInt *vals;
1366e26ddf31SBarry Smith 
1367e26ddf31SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_INT_ARRAY,&amsopt);CHKERRQ(ierr);
13686e02237eSPeter Brune     ierr = PetscMalloc1((*n),(PetscInt**)&amsopt->data);CHKERRQ(ierr);
1369e26ddf31SBarry Smith     vals = (PetscInt*)amsopt->data;
1370e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1371e26ddf31SBarry Smith     amsopt->arraylength = *n;
1372e26ddf31SBarry Smith   }
137353acd3b1SBarry Smith   ierr = PetscOptionsGetIntArray(PetscOptionsObject.prefix,opt,value,n,set);CHKERRQ(ierr);
137461b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
137553acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%d",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,value[0]);CHKERRQ(ierr);
137653acd3b1SBarry Smith     for (i=1; i<*n; i++) {
137753acd3b1SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,",%d",value[i]);CHKERRQ(ierr);
137853acd3b1SBarry Smith     }
13792aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
138053acd3b1SBarry Smith   }
138153acd3b1SBarry Smith   PetscFunctionReturn(0);
138253acd3b1SBarry Smith }
138353acd3b1SBarry Smith 
138453acd3b1SBarry Smith #undef __FUNCT__
138553acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsStringArray"
138653acd3b1SBarry Smith /*@C
138753acd3b1SBarry Smith    PetscOptionsStringArray - Gets an array of string values for a particular
138853acd3b1SBarry Smith    option in the database. The values must be separated with commas with
138953acd3b1SBarry Smith    no intervening spaces.
139053acd3b1SBarry Smith 
13913f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
139253acd3b1SBarry Smith 
139353acd3b1SBarry Smith    Input Parameters:
139453acd3b1SBarry Smith +  opt - the option one is seeking
139553acd3b1SBarry Smith .  text - short string describing option
139653acd3b1SBarry Smith .  man - manual page for option
139753acd3b1SBarry Smith -  nmax - maximum number of strings
139853acd3b1SBarry Smith 
139953acd3b1SBarry Smith    Output Parameter:
140053acd3b1SBarry Smith +  value - location to copy strings
140153acd3b1SBarry Smith .  nmax - actual number of strings found
140253acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
140353acd3b1SBarry Smith 
140453acd3b1SBarry Smith    Level: beginner
140553acd3b1SBarry Smith 
140653acd3b1SBarry Smith    Notes:
140753acd3b1SBarry Smith    The user should pass in an array of pointers to char, to hold all the
140853acd3b1SBarry Smith    strings returned by this function.
140953acd3b1SBarry Smith 
141053acd3b1SBarry Smith    The user is responsible for deallocating the strings that are
141153acd3b1SBarry Smith    returned. The Fortran interface for this routine is not supported.
141253acd3b1SBarry Smith 
141353acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
141453acd3b1SBarry Smith 
141553acd3b1SBarry Smith    Concepts: options database^array of strings
141653acd3b1SBarry Smith 
141753acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1418acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
141953acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
142053acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1421acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1422a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
142353acd3b1SBarry Smith @*/
14247087cfbeSBarry Smith PetscErrorCode  PetscOptionsStringArray(const char opt[],const char text[],const char man[],char *value[],PetscInt *nmax,PetscBool  *set)
142553acd3b1SBarry Smith {
142653acd3b1SBarry Smith   PetscErrorCode ierr;
14271ae3d29cSBarry Smith   PetscOptions   amsopt;
142853acd3b1SBarry Smith 
142953acd3b1SBarry Smith   PetscFunctionBegin;
14301ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
14311ae3d29cSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_STRING_ARRAY,&amsopt);CHKERRQ(ierr);
14326e02237eSPeter Brune     ierr = PetscMalloc1((*nmax),(char**)&amsopt->data);CHKERRQ(ierr);
1433a297a907SKarl Rupp 
14341ae3d29cSBarry Smith     amsopt->arraylength = *nmax;
14351ae3d29cSBarry Smith   }
143653acd3b1SBarry Smith   ierr = PetscOptionsGetStringArray(PetscOptionsObject.prefix,opt,value,nmax,set);CHKERRQ(ierr);
143761b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
14382aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <string1,string2,...>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
143953acd3b1SBarry Smith   }
144053acd3b1SBarry Smith   PetscFunctionReturn(0);
144153acd3b1SBarry Smith }
144253acd3b1SBarry Smith 
1443e2446a98SMatthew Knepley #undef __FUNCT__
1444acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolArray"
1445e2446a98SMatthew Knepley /*@C
1446acfcf0e5SJed Brown    PetscOptionsBoolArray - Gets an array of logical values (true or false) for a particular
1447e2446a98SMatthew Knepley    option in the database. The values must be separated with commas with
1448e2446a98SMatthew Knepley    no intervening spaces.
1449e2446a98SMatthew Knepley 
14503f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
1451e2446a98SMatthew Knepley 
1452e2446a98SMatthew Knepley    Input Parameters:
1453e2446a98SMatthew Knepley +  opt - the option one is seeking
1454e2446a98SMatthew Knepley .  text - short string describing option
1455e2446a98SMatthew Knepley .  man - manual page for option
1456e2446a98SMatthew Knepley -  nmax - maximum number of values
1457e2446a98SMatthew Knepley 
1458e2446a98SMatthew Knepley    Output Parameter:
1459e2446a98SMatthew Knepley +  value - location to copy values
1460e2446a98SMatthew Knepley .  nmax - actual number of values found
1461e2446a98SMatthew Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
1462e2446a98SMatthew Knepley 
1463e2446a98SMatthew Knepley    Level: beginner
1464e2446a98SMatthew Knepley 
1465e2446a98SMatthew Knepley    Notes:
1466e2446a98SMatthew Knepley    The user should pass in an array of doubles
1467e2446a98SMatthew Knepley 
1468e2446a98SMatthew Knepley    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
1469e2446a98SMatthew Knepley 
1470e2446a98SMatthew Knepley    Concepts: options database^array of strings
1471e2446a98SMatthew Knepley 
1472e2446a98SMatthew Knepley .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1473acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
1474e2446a98SMatthew Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1475e2446a98SMatthew Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1476acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1477a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
1478e2446a98SMatthew Knepley @*/
14797087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolArray(const char opt[],const char text[],const char man[],PetscBool value[],PetscInt *n,PetscBool *set)
1480e2446a98SMatthew Knepley {
1481e2446a98SMatthew Knepley   PetscErrorCode ierr;
1482e2446a98SMatthew Knepley   PetscInt       i;
14831ae3d29cSBarry Smith   PetscOptions   amsopt;
1484e2446a98SMatthew Knepley 
1485e2446a98SMatthew Knepley   PetscFunctionBegin;
14861ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
1487ace3abfcSBarry Smith     PetscBool *vals;
14881ae3d29cSBarry Smith 
14897781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL_ARRAY,&amsopt);CHKERRQ(ierr);
14906e02237eSPeter Brune     ierr = PetscMalloc1((*n),(PetscBool**)&amsopt->data);CHKERRQ(ierr);
1491ace3abfcSBarry Smith     vals = (PetscBool*)amsopt->data;
14921ae3d29cSBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
14931ae3d29cSBarry Smith     amsopt->arraylength = *n;
14941ae3d29cSBarry Smith   }
1495acfcf0e5SJed Brown   ierr = PetscOptionsGetBoolArray(PetscOptionsObject.prefix,opt,value,n,set);CHKERRQ(ierr);
1496e2446a98SMatthew Knepley   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
1497e2446a98SMatthew Knepley     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%d",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,value[0]);CHKERRQ(ierr);
1498e2446a98SMatthew Knepley     for (i=1; i<*n; i++) {
1499e2446a98SMatthew Knepley       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,",%d",value[i]);CHKERRQ(ierr);
1500e2446a98SMatthew Knepley     }
15012aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
1502e2446a98SMatthew Knepley   }
1503e2446a98SMatthew Knepley   PetscFunctionReturn(0);
1504e2446a98SMatthew Knepley }
1505e2446a98SMatthew Knepley 
15068cc676e6SMatthew G Knepley #undef __FUNCT__
15078cc676e6SMatthew G Knepley #define __FUNCT__ "PetscOptionsViewer"
15088cc676e6SMatthew G Knepley /*@C
15098cc676e6SMatthew G Knepley    PetscOptionsInt - Gets a viewer appropriate for the type indicated by the user
15108cc676e6SMatthew G Knepley 
15118cc676e6SMatthew G Knepley    Logically Collective on the communicator passed in PetscOptionsBegin()
15128cc676e6SMatthew G Knepley 
15138cc676e6SMatthew G Knepley    Input Parameters:
15148cc676e6SMatthew G Knepley +  opt - option name
15158cc676e6SMatthew G Knepley .  text - short string that describes the option
15168cc676e6SMatthew G Knepley -  man - manual page with additional information on option
15178cc676e6SMatthew G Knepley 
15188cc676e6SMatthew G Knepley    Output Parameter:
15198cc676e6SMatthew G Knepley +  viewer - the viewer
15208cc676e6SMatthew G Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
15218cc676e6SMatthew G Knepley 
15228cc676e6SMatthew G Knepley    Level: beginner
15238cc676e6SMatthew G Knepley 
15248cc676e6SMatthew G Knepley    Concepts: options database^has int
15258cc676e6SMatthew G Knepley 
15268cc676e6SMatthew G Knepley    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
15278cc676e6SMatthew G Knepley If no value is provided ascii:stdout is used
15288cc676e6SMatthew G Knepley $       ascii[:[filename][:format]]   defaults to stdout - format can be one of info, info_detailed, or matlab, for example ascii::info prints just the info
15298cc676e6SMatthew G Knepley $                                     about the object to standard out
15308cc676e6SMatthew G Knepley $       binary[:filename]   defaults to binaryoutput
15318cc676e6SMatthew G Knepley $       draw
15328cc676e6SMatthew G Knepley $       socket[:port]    defaults to the standard output port
15338cc676e6SMatthew G Knepley 
1534cffb1e40SBarry Smith    Use PetscRestoreViewerDestroy() after using the viewer, otherwise a memory leak will occur
15358cc676e6SMatthew G Knepley 
15368cc676e6SMatthew G Knepley .seealso: PetscOptionsGetViewer(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
15378cc676e6SMatthew G Knepley           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
15388cc676e6SMatthew G Knepley           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
15398cc676e6SMatthew G Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
15408cc676e6SMatthew G Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
15418cc676e6SMatthew G Knepley           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1542a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
15438cc676e6SMatthew G Knepley @*/
1544cffb1e40SBarry Smith PetscErrorCode  PetscOptionsViewer(const char opt[],const char text[],const char man[],PetscViewer *viewer,PetscViewerFormat *format,PetscBool  *set)
15458cc676e6SMatthew G Knepley {
15468cc676e6SMatthew G Knepley   PetscErrorCode ierr;
15478cc676e6SMatthew G Knepley   PetscOptions   amsopt;
15488cc676e6SMatthew G Knepley 
15498cc676e6SMatthew G Knepley   PetscFunctionBegin;
15508cc676e6SMatthew G Knepley   if (!PetscOptionsPublishCount) {
15518cc676e6SMatthew G Knepley     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
155264facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
15535b02f95dSBarry Smith     ierr = PetscStrdup("",(char**)&amsopt->data);CHKERRQ(ierr);
15548cc676e6SMatthew G Knepley   }
1555cffb1e40SBarry Smith   ierr = PetscOptionsGetViewer(PetscOptionsObject.comm,PetscOptionsObject.prefix,opt,viewer,format,set);CHKERRQ(ierr);
15568cc676e6SMatthew G Knepley   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
15578cc676e6SMatthew G Knepley     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%s>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,"",text,ManSection(man));CHKERRQ(ierr);
15588cc676e6SMatthew G Knepley   }
15598cc676e6SMatthew G Knepley   PetscFunctionReturn(0);
15608cc676e6SMatthew G Knepley }
15618cc676e6SMatthew G Knepley 
156253acd3b1SBarry Smith 
156353acd3b1SBarry Smith #undef __FUNCT__
156453acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsHead"
156553acd3b1SBarry Smith /*@C
1566b52f573bSBarry Smith      PetscOptionsHead - Puts a heading before listing any more published options. Used, for example,
156753acd3b1SBarry Smith             in KSPSetFromOptions_GMRES().
156853acd3b1SBarry Smith 
15693f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
157053acd3b1SBarry Smith 
157153acd3b1SBarry Smith    Input Parameter:
157253acd3b1SBarry Smith .   head - the heading text
157353acd3b1SBarry Smith 
157453acd3b1SBarry Smith 
157553acd3b1SBarry Smith    Level: intermediate
157653acd3b1SBarry Smith 
157753acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
157853acd3b1SBarry Smith 
1579b52f573bSBarry Smith           Can be followed by a call to PetscOptionsTail() in the same function.
158053acd3b1SBarry Smith 
158153acd3b1SBarry Smith    Concepts: options database^subheading
158253acd3b1SBarry Smith 
158353acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1584acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
158553acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
158653acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1587acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1588a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
158953acd3b1SBarry Smith @*/
15907087cfbeSBarry Smith PetscErrorCode  PetscOptionsHead(const char head[])
159153acd3b1SBarry Smith {
159253acd3b1SBarry Smith   PetscErrorCode ierr;
159353acd3b1SBarry Smith 
159453acd3b1SBarry Smith   PetscFunctionBegin;
159561b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
159653acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  %s\n",head);CHKERRQ(ierr);
159753acd3b1SBarry Smith   }
159853acd3b1SBarry Smith   PetscFunctionReturn(0);
159953acd3b1SBarry Smith }
160053acd3b1SBarry Smith 
160153acd3b1SBarry Smith 
160253acd3b1SBarry Smith 
160353acd3b1SBarry Smith 
160453acd3b1SBarry Smith 
160553acd3b1SBarry Smith 
1606