xref: /petsc/src/sys/objects/aoptions.c (revision 9c1e0ce8770ab6d4884a7bf923ea63b737ed5291)
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*9c1e0ce8SBarry Smith /*
37764bbc9efSBarry Smith                                    "<script jowererw type=\"text/javascript\" src=\"https://bitbucket.org/saws/saws/raw/master/js/jquery-1.9.1.js\"></script>\n"
37864bbc9efSBarry Smith                                    "<script type=\"text/javascript\" src=\"https://bitbucket.org/saws/saws/raw/master/js/jsSAWs.js\"></script>\n"
379*9c1e0ce8SBarry Smith  */
380*9c1e0ce8SBarry Smith 
381*9c1e0ce8SBarry Smith static const char *OptionsHeader = "<head>\n"
382*9c1e0ce8SBarry Smith                                    "<script type=\"text/javascript\" src=\"js/jquery-1.9.1.js\"></script>\n"
383*9c1e0ce8SBarry Smith                                    "<script type=\"text/javascript\" src=\"js/jsSAWs.js\"></script>\n"
384d1fc0251SBarry Smith                                    "<script type=\"text/javascript\" src=\"js/PETSc.js\"></script>\n"
38564bbc9efSBarry Smith                                    "<script>\n"
38664bbc9efSBarry Smith                                       "jQuery(document).ready(function() {\n"
38764bbc9efSBarry Smith                                          "PETSc.getAndDisplayDirectory(null,\"#variablesInfo\")\n"
38864bbc9efSBarry Smith                                       "})\n"
38964bbc9efSBarry Smith                                   "</script>\n"
39064bbc9efSBarry Smith                                   "</head>\n";
39164bbc9efSBarry Smith static const char *OptionsBodyBottom = "<div id=\"variablesInfo\" style=\"float:left\"></div>\n<br>\n</body>";
39264bbc9efSBarry Smith 
393d5649816SBarry Smith #undef __FUNCT__
3947781c08eSBarry Smith #define __FUNCT__ "PetscOptionsSAWsInput"
395b3506946SBarry Smith /*
3967781c08eSBarry Smith     PetscOptionsSAWsInput - Presents all the PETSc Options processed by the program so the user may change them at runtime using the SAWs
397b3506946SBarry Smith 
398b3506946SBarry Smith     Bugs:
399b3506946SBarry Smith +    All processes must traverse through the exact same set of option queries do to the call to PetscScanString()
400b3506946SBarry Smith .    Internal strings have arbitrary length and string copies are not checked that they fit into string space
401b3506946SBarry Smith -    Only works for PetscInt == int, PetscReal == double etc
402b3506946SBarry Smith 
403b3506946SBarry Smith 
404b3506946SBarry Smith */
405475446a1SBarry Smith PetscErrorCode PetscOptionsSAWsInput()
406b3506946SBarry Smith {
407b3506946SBarry Smith   PetscErrorCode ierr;
408b3506946SBarry Smith   PetscOptions   next     = PetscOptionsObject.next;
409d5649816SBarry Smith   static int     mancount = 0;
410b3506946SBarry Smith   char           options[16];
4117aab2a10SBarry Smith   PetscBool      changedmethod = PETSC_FALSE;
41288a9752cSBarry Smith   char           manname[16],textname[16];
4132657e9d9SBarry Smith   char           dir[1024];
414b3506946SBarry Smith 
415b3506946SBarry Smith   /* the next line is a bug, this will only work if all processors are here, the comm passed in is ignored!!! */
416b3506946SBarry Smith   sprintf(options,"Options_%d",count++);
417a297a907SKarl Rupp 
418e04113cfSBarry Smith   PetscOptionsObject.pprefix = PetscOptionsObject.prefix; /* SAWs will change this, so cannot pass prefix directly */
4191bc75a8dSBarry Smith 
420d50831c4SBarry Smith   ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s","_title");CHKERRQ(ierr);
4217781c08eSBarry Smith   PetscStackCallSAWs(SAWs_Register,(dir,&PetscOptionsObject.title,1,SAWs_READ,SAWs_STRING));
4227781c08eSBarry Smith   ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s","prefix");CHKERRQ(ierr);
4232657e9d9SBarry Smith   PetscStackCallSAWs(SAWs_Register,(dir,&PetscOptionsObject.pprefix,1,SAWs_READ,SAWs_STRING));
4242657e9d9SBarry Smith   PetscStackCallSAWs(SAWs_Register,("/PETSc/Options/ChangedMethod",&changedmethod,1,SAWs_WRITE,SAWs_BOOLEAN));
425b3506946SBarry Smith 
426b3506946SBarry Smith   while (next) {
427d50831c4SBarry Smith     sprintf(manname,"_man_%d",mancount);
4282657e9d9SBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",manname);CHKERRQ(ierr);
4297781c08eSBarry Smith     PetscStackCallSAWs(SAWs_Register,(dir,&next->man,1,SAWs_READ,SAWs_STRING));
430d50831c4SBarry Smith     sprintf(textname,"_text_%d",mancount++);
4317781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",textname);CHKERRQ(ierr);
4327781c08eSBarry Smith     PetscStackCallSAWs(SAWs_Register,(dir,&next->text,1,SAWs_READ,SAWs_STRING));
4339f32e415SBarry Smith 
434b3506946SBarry Smith     switch (next->type) {
435b3506946SBarry Smith     case OPTION_HEAD:
436b3506946SBarry Smith       break;
437b3506946SBarry Smith     case OPTION_INT_ARRAY:
4387781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4392657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_INT));
440b3506946SBarry Smith       break;
441b3506946SBarry Smith     case OPTION_REAL_ARRAY:
4427781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4432657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_DOUBLE));
444b3506946SBarry Smith       break;
445b3506946SBarry Smith     case OPTION_INT:
4467781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4472657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_INT));
448b3506946SBarry Smith       break;
449b3506946SBarry Smith     case OPTION_REAL:
4507781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4512657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_DOUBLE));
452b3506946SBarry Smith       break;
4537781c08eSBarry Smith     case OPTION_BOOL:
4547781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4552657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_BOOLEAN));
4561ae3d29cSBarry Smith       break;
4577781c08eSBarry Smith     case OPTION_BOOL_ARRAY:
4587781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4592657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_BOOLEAN));
46071f08665SBarry Smith       break;
461b3506946SBarry Smith     case OPTION_STRING:
4627781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4637781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
4641ae3d29cSBarry Smith       break;
4651ae3d29cSBarry Smith     case OPTION_STRING_ARRAY:
4667781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4672657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_STRING));
468b3506946SBarry Smith       break;
469a264d7a6SBarry Smith     case OPTION_FLIST:
470a264d7a6SBarry Smith       {
471a264d7a6SBarry Smith       PetscInt ntext;
4727781c08eSBarry Smith       ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4737781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
474a264d7a6SBarry Smith       ierr = PetscFunctionListGet(next->flist,(const char***)&next->edata,&ntext);CHKERRQ(ierr);
475a264d7a6SBarry Smith       PetscStackCallSAWs(SAWs_Set_Legal_Variable_Values,(dir,ntext,next->edata));
476a264d7a6SBarry Smith       }
477a264d7a6SBarry Smith       break;
4781ae3d29cSBarry Smith     case OPTION_ELIST:
479a264d7a6SBarry Smith       {
480a264d7a6SBarry Smith       PetscInt ntext = next->nlist;
4817781c08eSBarry Smith       ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4827781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
483ead66b60SBarry Smith       ierr = PetscMalloc1((ntext+1),(char***)&next->edata);CHKERRQ(ierr);
4841ae3d29cSBarry Smith       ierr = PetscMemcpy(next->edata,next->list,ntext*sizeof(char*));CHKERRQ(ierr);
485a264d7a6SBarry Smith       PetscStackCallSAWs(SAWs_Set_Legal_Variable_Values,(dir,ntext,next->edata));
486a264d7a6SBarry Smith       }
487a264d7a6SBarry Smith       break;
488b3506946SBarry Smith     default:
489b3506946SBarry Smith       break;
490b3506946SBarry Smith     }
491b3506946SBarry Smith     next = next->next;
492b3506946SBarry Smith   }
493b3506946SBarry Smith 
494b3506946SBarry Smith   /* wait until accessor has unlocked the memory */
49564bbc9efSBarry Smith   PetscStackCallSAWs(SAWs_Push_Header,("index.html",OptionsHeader));
49664bbc9efSBarry Smith   PetscStackCallSAWs(SAWs_Push_Body,("index.html",2,OptionsBodyBottom));
4977aab2a10SBarry Smith   ierr = PetscSAWsBlock();CHKERRQ(ierr);
49864bbc9efSBarry Smith   PetscStackCallSAWs(SAWs_Pop_Header,("index.html"));
49964bbc9efSBarry Smith   PetscStackCallSAWs(SAWs_Pop_Body,("index.html",2));
500b3506946SBarry Smith 
50188a9752cSBarry Smith   /* determine if any values have been set in GUI */
50288a9752cSBarry Smith   next = PetscOptionsObject.next;
50388a9752cSBarry Smith   while (next) {
50488a9752cSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
50588a9752cSBarry Smith     PetscStackCallSAWs(SAWs_Selected,(dir,&next->set));
50688a9752cSBarry Smith     next = next->next;
50788a9752cSBarry Smith   }
50888a9752cSBarry Smith 
509b3506946SBarry Smith   /* reset counter to -2; this updates the screen with the new options for the selected method */
510b3506946SBarry Smith   if (changedmethod) PetscOptionsPublishCount = -2;
511b3506946SBarry Smith 
5129a492a5cSBarry Smith   PetscStackCallSAWs(SAWs_Delete,("/PETSc/Options"));
513b3506946SBarry Smith   PetscFunctionReturn(0);
514b3506946SBarry Smith }
515b3506946SBarry Smith #endif
516b3506946SBarry Smith 
5176356e834SBarry Smith #undef __FUNCT__
51853acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEnd_Private"
51953acd3b1SBarry Smith PetscErrorCode PetscOptionsEnd_Private(void)
52053acd3b1SBarry Smith {
52153acd3b1SBarry Smith   PetscErrorCode ierr;
5226356e834SBarry Smith   PetscOptions   last;
5236356e834SBarry Smith   char           option[256],value[1024],tmp[32];
524330cf3c9SBarry Smith   size_t         j;
52553acd3b1SBarry Smith 
52653acd3b1SBarry Smith   PetscFunctionBegin;
527aee2cecaSBarry Smith   if (PetscOptionsObject.next) {
528b3506946SBarry Smith     if (!PetscOptionsPublishCount) {
529a264d7a6SBarry Smith #if defined(PETSC_HAVE_SAWS)
530475446a1SBarry Smith       ierr = PetscOptionsSAWsInput();CHKERRQ(ierr);
531b3506946SBarry Smith #else
53271f08665SBarry Smith       ierr = PetscOptionsGetFromTextInput();CHKERRQ(ierr);
533b3506946SBarry Smith #endif
534aee2cecaSBarry Smith     }
535aee2cecaSBarry Smith   }
5366356e834SBarry Smith 
537c31cb41cSBarry Smith   ierr = PetscFree(PetscOptionsObject.title);CHKERRQ(ierr);
538c31cb41cSBarry Smith   ierr = PetscFree(PetscOptionsObject.prefix);CHKERRQ(ierr);
5396356e834SBarry Smith 
5406356e834SBarry Smith   /* reset counter to -2; this updates the screen with the new options for the selected method */
5416356e834SBarry Smith   if (PetscOptionsObject.changedmethod) PetscOptionsPublishCount = -2;
54261b37b28SSatish Balay   /* reset alreadyprinted flag */
54361b37b28SSatish Balay   PetscOptionsObject.alreadyprinted = PETSC_FALSE;
5443194b578SJed Brown   if (PetscOptionsObject.object) PetscOptionsObject.object->optionsprinted = PETSC_TRUE;
5450298fd71SBarry Smith   PetscOptionsObject.object = NULL;
5466356e834SBarry Smith 
5476356e834SBarry Smith   while (PetscOptionsObject.next) {
5486356e834SBarry Smith     if (PetscOptionsObject.next->set) {
5496356e834SBarry Smith       if (PetscOptionsObject.prefix) {
5506356e834SBarry Smith         ierr = PetscStrcpy(option,"-");CHKERRQ(ierr);
5516356e834SBarry Smith         ierr = PetscStrcat(option,PetscOptionsObject.prefix);CHKERRQ(ierr);
5526356e834SBarry Smith         ierr = PetscStrcat(option,PetscOptionsObject.next->option+1);CHKERRQ(ierr);
5536356e834SBarry Smith       } else {
5546356e834SBarry Smith         ierr = PetscStrcpy(option,PetscOptionsObject.next->option);CHKERRQ(ierr);
5556356e834SBarry Smith       }
5566356e834SBarry Smith 
5576356e834SBarry Smith       switch (PetscOptionsObject.next->type) {
5586356e834SBarry Smith       case OPTION_HEAD:
5596356e834SBarry Smith         break;
560e26ddf31SBarry Smith       case OPTION_INT_ARRAY:
561e26ddf31SBarry Smith         sprintf(value,"%d",(int)((PetscInt*)PetscOptionsObject.next->data)[0]);
562e26ddf31SBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
563e26ddf31SBarry Smith           sprintf(tmp,"%d",(int)((PetscInt*)PetscOptionsObject.next->data)[j]);
564e26ddf31SBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
565e26ddf31SBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
566e26ddf31SBarry Smith         }
567e26ddf31SBarry Smith         break;
5686356e834SBarry Smith       case OPTION_INT:
5697a72a596SBarry Smith         sprintf(value,"%d",(int) *(PetscInt*)PetscOptionsObject.next->data);
5706356e834SBarry Smith         break;
5716356e834SBarry Smith       case OPTION_REAL:
5727a72a596SBarry Smith         sprintf(value,"%g",(double) *(PetscReal*)PetscOptionsObject.next->data);
5736356e834SBarry Smith         break;
5746356e834SBarry Smith       case OPTION_REAL_ARRAY:
5757a72a596SBarry Smith         sprintf(value,"%g",(double)((PetscReal*)PetscOptionsObject.next->data)[0]);
5766356e834SBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
5777a72a596SBarry Smith           sprintf(tmp,"%g",(double)((PetscReal*)PetscOptionsObject.next->data)[j]);
5786356e834SBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5796356e834SBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5806356e834SBarry Smith         }
5816356e834SBarry Smith         break;
5827781c08eSBarry Smith       case OPTION_BOOL:
58371f08665SBarry Smith         sprintf(value,"%d",*(int*)PetscOptionsObject.next->data);
5846356e834SBarry Smith         break;
5857781c08eSBarry Smith       case OPTION_BOOL_ARRAY:
586ace3abfcSBarry Smith         sprintf(value,"%d",(int)((PetscBool*)PetscOptionsObject.next->data)[0]);
5871ae3d29cSBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
588ace3abfcSBarry Smith           sprintf(tmp,"%d",(int)((PetscBool*)PetscOptionsObject.next->data)[j]);
5891ae3d29cSBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5901ae3d29cSBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5911ae3d29cSBarry Smith         }
5921ae3d29cSBarry Smith         break;
593a264d7a6SBarry Smith       case OPTION_FLIST:
5941ae3d29cSBarry Smith       case OPTION_ELIST:
5957781c08eSBarry Smith         ierr = PetscStrcpy(value,(char*)PetscOptionsObject.next->data);CHKERRQ(ierr);
5966356e834SBarry Smith         break;
5971ae3d29cSBarry Smith       case OPTION_STRING:
598475446a1SBarry Smith         ierr = PetscStrcpy(value,(char*)PetscOptionsObject.next->data);CHKERRQ(ierr);
5998da2146eSBarry Smith         break;
6001ae3d29cSBarry Smith       case OPTION_STRING_ARRAY:
6011ae3d29cSBarry Smith         sprintf(value,"%s",((char**)PetscOptionsObject.next->data)[0]);
6021ae3d29cSBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
6031ae3d29cSBarry Smith           sprintf(tmp,"%s",((char**)PetscOptionsObject.next->data)[j]);
6041ae3d29cSBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
6051ae3d29cSBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
6061ae3d29cSBarry Smith         }
6076356e834SBarry Smith         break;
6086356e834SBarry Smith       }
6096356e834SBarry Smith       ierr = PetscOptionsSetValue(option,value);CHKERRQ(ierr);
6106356e834SBarry Smith     }
611503cfb0cSBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->text);CHKERRQ(ierr);
612503cfb0cSBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->option);CHKERRQ(ierr);
6136356e834SBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->man);CHKERRQ(ierr);
61471f08665SBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->edata);CHKERRQ(ierr);
615c979a496SBarry Smith 
616c979a496SBarry Smith     if ((PetscOptionsObject.next->type == OPTION_STRING) || (PetscOptionsObject.next->type == OPTION_FLIST) || (PetscOptionsObject.next->type == OPTION_ELIST)){
61764facd6cSBarry Smith       /* must use system free since SAWs may have allocated it */
618c979a496SBarry Smith       free(PetscOptionsObject.next->data);
619c979a496SBarry Smith     } else {
6207781c08eSBarry Smith       ierr   = PetscFree(PetscOptionsObject.next->data);CHKERRQ(ierr);
621c979a496SBarry Smith     }
6227781c08eSBarry Smith 
6236356e834SBarry Smith     last                    = PetscOptionsObject.next;
6246356e834SBarry Smith     PetscOptionsObject.next = PetscOptionsObject.next->next;
6256356e834SBarry Smith     ierr                    = PetscFree(last);CHKERRQ(ierr);
6266356e834SBarry Smith   }
6276356e834SBarry Smith   PetscOptionsObject.next = 0;
62853acd3b1SBarry Smith   PetscFunctionReturn(0);
62953acd3b1SBarry Smith }
63053acd3b1SBarry Smith 
63153acd3b1SBarry Smith #undef __FUNCT__
63253acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEnum"
63353acd3b1SBarry Smith /*@C
63453acd3b1SBarry Smith    PetscOptionsEnum - Gets the enum value for a particular option in the database.
63553acd3b1SBarry Smith 
6363f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
63753acd3b1SBarry Smith 
63853acd3b1SBarry Smith    Input Parameters:
63953acd3b1SBarry Smith +  opt - option name
64053acd3b1SBarry Smith .  text - short string that describes the option
64153acd3b1SBarry Smith .  man - manual page with additional information on option
64253acd3b1SBarry Smith .  list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null
64353acd3b1SBarry Smith -  defaultv - the default (current) value
64453acd3b1SBarry Smith 
64553acd3b1SBarry Smith    Output Parameter:
64653acd3b1SBarry Smith +  value - the  value to return
647b32e0204SMatthew G Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
64853acd3b1SBarry Smith 
64953acd3b1SBarry Smith    Level: beginner
65053acd3b1SBarry Smith 
65153acd3b1SBarry Smith    Concepts: options database
65253acd3b1SBarry Smith 
65353acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
65453acd3b1SBarry Smith 
65553acd3b1SBarry Smith           list is usually something like PCASMTypes or some other predefined list of enum names
65653acd3b1SBarry Smith 
65753acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
658acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
659acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
66053acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
66153acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
662acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
663a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
66453acd3b1SBarry Smith @*/
6657087cfbeSBarry Smith PetscErrorCode  PetscOptionsEnum(const char opt[],const char text[],const char man[],const char *const *list,PetscEnum defaultv,PetscEnum *value,PetscBool  *set)
66653acd3b1SBarry Smith {
66753acd3b1SBarry Smith   PetscErrorCode ierr;
66853acd3b1SBarry Smith   PetscInt       ntext = 0;
669aa5bb8c0SSatish Balay   PetscInt       tval;
670ace3abfcSBarry Smith   PetscBool      tflg;
67153acd3b1SBarry Smith 
67253acd3b1SBarry Smith   PetscFunctionBegin;
67353acd3b1SBarry Smith   while (list[ntext++]) {
674e32f2f54SBarry Smith     if (ntext > 50) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument appears to be wrong or have more than 50 entries");
67553acd3b1SBarry Smith   }
676e32f2f54SBarry Smith   if (ntext < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument must have at least two entries: typename and type prefix");
67753acd3b1SBarry Smith   ntext -= 3;
678aa5bb8c0SSatish Balay   ierr   = PetscOptionsEList(opt,text,man,list,ntext,list[defaultv],&tval,&tflg);CHKERRQ(ierr);
679aa5bb8c0SSatish Balay   /* with PETSC_USE_64BIT_INDICES sizeof(PetscInt) != sizeof(PetscEnum) */
680aa5bb8c0SSatish Balay   if (tflg) *value = (PetscEnum)tval;
681aa5bb8c0SSatish Balay   if (set)  *set   = tflg;
68253acd3b1SBarry Smith   PetscFunctionReturn(0);
68353acd3b1SBarry Smith }
68453acd3b1SBarry Smith 
68553acd3b1SBarry Smith /* -------------------------------------------------------------------------------------------------------------*/
68653acd3b1SBarry Smith #undef __FUNCT__
68753acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsInt"
68853acd3b1SBarry Smith /*@C
68953acd3b1SBarry Smith    PetscOptionsInt - Gets the integer value for a particular option in the database.
69053acd3b1SBarry Smith 
6913f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
69253acd3b1SBarry Smith 
69353acd3b1SBarry Smith    Input Parameters:
69453acd3b1SBarry Smith +  opt - option name
69553acd3b1SBarry Smith .  text - short string that describes the option
69653acd3b1SBarry Smith .  man - manual page with additional information on option
69753acd3b1SBarry Smith -  defaultv - the default (current) value
69853acd3b1SBarry Smith 
69953acd3b1SBarry Smith    Output Parameter:
70053acd3b1SBarry Smith +  value - the integer value to return
70153acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
70253acd3b1SBarry Smith 
70353acd3b1SBarry Smith    Level: beginner
70453acd3b1SBarry Smith 
70553acd3b1SBarry Smith    Concepts: options database^has int
70653acd3b1SBarry Smith 
70753acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
70853acd3b1SBarry Smith 
70953acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
710acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
711acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
71253acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
71353acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
714acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
715a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
71653acd3b1SBarry Smith @*/
7177087cfbeSBarry Smith PetscErrorCode  PetscOptionsInt(const char opt[],const char text[],const char man[],PetscInt defaultv,PetscInt *value,PetscBool  *set)
71853acd3b1SBarry Smith {
71953acd3b1SBarry Smith   PetscErrorCode ierr;
7206356e834SBarry Smith   PetscOptions   amsopt;
72153acd3b1SBarry Smith 
72253acd3b1SBarry Smith   PetscFunctionBegin;
723b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
7246356e834SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_INT,&amsopt);CHKERRQ(ierr);
7256356e834SBarry Smith     ierr = PetscMalloc(sizeof(PetscInt),&amsopt->data);CHKERRQ(ierr);
726a297a907SKarl Rupp 
7276356e834SBarry Smith     *(PetscInt*)amsopt->data = defaultv;
728af6d86caSBarry Smith   }
72953acd3b1SBarry Smith   ierr = PetscOptionsGetInt(PetscOptionsObject.prefix,opt,value,set);CHKERRQ(ierr);
73061b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
7312aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%d>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,defaultv,text,ManSection(man));CHKERRQ(ierr);
73253acd3b1SBarry Smith   }
73353acd3b1SBarry Smith   PetscFunctionReturn(0);
73453acd3b1SBarry Smith }
73553acd3b1SBarry Smith 
73653acd3b1SBarry Smith #undef __FUNCT__
73753acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsString"
73853acd3b1SBarry Smith /*@C
73953acd3b1SBarry Smith    PetscOptionsString - Gets the string value for a particular option in the database.
74053acd3b1SBarry Smith 
7413f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
74253acd3b1SBarry Smith 
74353acd3b1SBarry Smith    Input Parameters:
74453acd3b1SBarry Smith +  opt - option name
74553acd3b1SBarry Smith .  text - short string that describes the option
74653acd3b1SBarry Smith .  man - manual page with additional information on option
747bcbf2dc5SJed Brown .  defaultv - the default (current) value
748bcbf2dc5SJed Brown -  len - length of the result string including null terminator
74953acd3b1SBarry Smith 
75053acd3b1SBarry Smith    Output Parameter:
75153acd3b1SBarry Smith +  value - the value to return
75253acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
75353acd3b1SBarry Smith 
75453acd3b1SBarry Smith    Level: beginner
75553acd3b1SBarry Smith 
75653acd3b1SBarry Smith    Concepts: options database^has int
75753acd3b1SBarry Smith 
75853acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
75953acd3b1SBarry Smith 
7607fccdfe4SBarry 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).
7617fccdfe4SBarry Smith 
76253acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
763acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
764acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsReal(), PetscOptionsBool(),
76553acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
76653acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
767acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
768a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
76953acd3b1SBarry Smith @*/
7707087cfbeSBarry Smith PetscErrorCode  PetscOptionsString(const char opt[],const char text[],const char man[],const char defaultv[],char value[],size_t len,PetscBool  *set)
77153acd3b1SBarry Smith {
77253acd3b1SBarry Smith   PetscErrorCode ierr;
773aee2cecaSBarry Smith   PetscOptions   amsopt;
77453acd3b1SBarry Smith 
77553acd3b1SBarry Smith   PetscFunctionBegin;
776b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
777aee2cecaSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
77864facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
7795b02f95dSBarry Smith     ierr = PetscStrdup(defaultv ? defaultv : "",(char**)&amsopt->data);CHKERRQ(ierr);
780af6d86caSBarry Smith   }
78153acd3b1SBarry Smith   ierr = PetscOptionsGetString(PetscOptionsObject.prefix,opt,value,len,set);CHKERRQ(ierr);
78261b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
7832aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%s>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,defaultv,text,ManSection(man));CHKERRQ(ierr);
78453acd3b1SBarry Smith   }
78553acd3b1SBarry Smith   PetscFunctionReturn(0);
78653acd3b1SBarry Smith }
78753acd3b1SBarry Smith 
78853acd3b1SBarry Smith #undef __FUNCT__
78953acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsReal"
79053acd3b1SBarry Smith /*@C
79153acd3b1SBarry Smith    PetscOptionsReal - Gets the PetscReal value for a particular option in the database.
79253acd3b1SBarry Smith 
7933f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
79453acd3b1SBarry Smith 
79553acd3b1SBarry Smith    Input Parameters:
79653acd3b1SBarry Smith +  opt - option name
79753acd3b1SBarry Smith .  text - short string that describes the option
79853acd3b1SBarry Smith .  man - manual page with additional information on option
79953acd3b1SBarry Smith -  defaultv - the default (current) value
80053acd3b1SBarry Smith 
80153acd3b1SBarry Smith    Output Parameter:
80253acd3b1SBarry Smith +  value - the value to return
80353acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
80453acd3b1SBarry Smith 
80553acd3b1SBarry Smith    Level: beginner
80653acd3b1SBarry Smith 
80753acd3b1SBarry Smith    Concepts: options database^has int
80853acd3b1SBarry Smith 
80953acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
81053acd3b1SBarry Smith 
81153acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
812acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
813acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
81453acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
81553acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
816acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
817a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
81853acd3b1SBarry Smith @*/
8197087cfbeSBarry Smith PetscErrorCode  PetscOptionsReal(const char opt[],const char text[],const char man[],PetscReal defaultv,PetscReal *value,PetscBool  *set)
82053acd3b1SBarry Smith {
82153acd3b1SBarry Smith   PetscErrorCode ierr;
822538aa990SBarry Smith   PetscOptions   amsopt;
82353acd3b1SBarry Smith 
82453acd3b1SBarry Smith   PetscFunctionBegin;
825b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
826538aa990SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_REAL,&amsopt);CHKERRQ(ierr);
827538aa990SBarry Smith     ierr = PetscMalloc(sizeof(PetscReal),&amsopt->data);CHKERRQ(ierr);
828a297a907SKarl Rupp 
829538aa990SBarry Smith     *(PetscReal*)amsopt->data = defaultv;
830538aa990SBarry Smith   }
83153acd3b1SBarry Smith   ierr = PetscOptionsGetReal(PetscOptionsObject.prefix,opt,value,set);CHKERRQ(ierr);
83261b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
8332aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%G>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,defaultv,text,ManSection(man));CHKERRQ(ierr);
83453acd3b1SBarry Smith   }
83553acd3b1SBarry Smith   PetscFunctionReturn(0);
83653acd3b1SBarry Smith }
83753acd3b1SBarry Smith 
83853acd3b1SBarry Smith #undef __FUNCT__
83953acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsScalar"
84053acd3b1SBarry Smith /*@C
84153acd3b1SBarry Smith    PetscOptionsScalar - Gets the scalar value for a particular option in the database.
84253acd3b1SBarry Smith 
8433f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
84453acd3b1SBarry Smith 
84553acd3b1SBarry Smith    Input Parameters:
84653acd3b1SBarry Smith +  opt - option name
84753acd3b1SBarry Smith .  text - short string that describes the option
84853acd3b1SBarry Smith .  man - manual page with additional information on option
84953acd3b1SBarry Smith -  defaultv - the default (current) value
85053acd3b1SBarry Smith 
85153acd3b1SBarry Smith    Output Parameter:
85253acd3b1SBarry Smith +  value - the value to return
85353acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
85453acd3b1SBarry Smith 
85553acd3b1SBarry Smith    Level: beginner
85653acd3b1SBarry Smith 
85753acd3b1SBarry Smith    Concepts: options database^has int
85853acd3b1SBarry Smith 
85953acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
86053acd3b1SBarry Smith 
86153acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
862acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
863acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
86453acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
86553acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
866acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
867a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
86853acd3b1SBarry Smith @*/
8697087cfbeSBarry Smith PetscErrorCode  PetscOptionsScalar(const char opt[],const char text[],const char man[],PetscScalar defaultv,PetscScalar *value,PetscBool  *set)
87053acd3b1SBarry Smith {
87153acd3b1SBarry Smith   PetscErrorCode ierr;
87253acd3b1SBarry Smith 
87353acd3b1SBarry Smith   PetscFunctionBegin;
87453acd3b1SBarry Smith #if !defined(PETSC_USE_COMPLEX)
87553acd3b1SBarry Smith   ierr = PetscOptionsReal(opt,text,man,defaultv,value,set);CHKERRQ(ierr);
87653acd3b1SBarry Smith #else
87753acd3b1SBarry Smith   ierr = PetscOptionsGetScalar(PetscOptionsObject.prefix,opt,value,set);CHKERRQ(ierr);
87853acd3b1SBarry Smith #endif
87953acd3b1SBarry Smith   PetscFunctionReturn(0);
88053acd3b1SBarry Smith }
88153acd3b1SBarry Smith 
88253acd3b1SBarry Smith #undef __FUNCT__
88353acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsName"
88453acd3b1SBarry Smith /*@C
88590d69ab7SBarry 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
88690d69ab7SBarry Smith                       its value is set to false.
88753acd3b1SBarry Smith 
8883f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
88953acd3b1SBarry Smith 
89053acd3b1SBarry Smith    Input Parameters:
89153acd3b1SBarry Smith +  opt - option name
89253acd3b1SBarry Smith .  text - short string that describes the option
89353acd3b1SBarry Smith -  man - manual page with additional information on option
89453acd3b1SBarry Smith 
89553acd3b1SBarry Smith    Output Parameter:
89653acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
89753acd3b1SBarry Smith 
89853acd3b1SBarry Smith    Level: beginner
89953acd3b1SBarry Smith 
90053acd3b1SBarry Smith    Concepts: options database^has int
90153acd3b1SBarry Smith 
90253acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
90353acd3b1SBarry Smith 
90453acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
905acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
906acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
90753acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
90853acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
909acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
910a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
91153acd3b1SBarry Smith @*/
9127087cfbeSBarry Smith PetscErrorCode  PetscOptionsName(const char opt[],const char text[],const char man[],PetscBool  *flg)
91353acd3b1SBarry Smith {
91453acd3b1SBarry Smith   PetscErrorCode ierr;
9151ae3d29cSBarry Smith   PetscOptions   amsopt;
91653acd3b1SBarry Smith 
91753acd3b1SBarry Smith   PetscFunctionBegin;
9181ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
9197781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
920ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
921a297a907SKarl Rupp 
922ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
9231ae3d29cSBarry Smith   }
92453acd3b1SBarry Smith   ierr = PetscOptionsHasName(PetscOptionsObject.prefix,opt,flg);CHKERRQ(ierr);
92561b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
9262aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
92753acd3b1SBarry Smith   }
92853acd3b1SBarry Smith   PetscFunctionReturn(0);
92953acd3b1SBarry Smith }
93053acd3b1SBarry Smith 
93153acd3b1SBarry Smith #undef __FUNCT__
932a264d7a6SBarry Smith #define __FUNCT__ "PetscOptionsFList"
93353acd3b1SBarry Smith /*@C
934a264d7a6SBarry Smith      PetscOptionsFList - Puts a list of option values that a single one may be selected from
93553acd3b1SBarry Smith 
9363f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
93753acd3b1SBarry Smith 
93853acd3b1SBarry Smith    Input Parameters:
93953acd3b1SBarry Smith +  opt - option name
94053acd3b1SBarry Smith .  text - short string that describes the option
94153acd3b1SBarry Smith .  man - manual page with additional information on option
94253acd3b1SBarry Smith .  list - the possible choices
9433cc1e11dSBarry Smith .  defaultv - the default (current) value
9443cc1e11dSBarry Smith -  len - the length of the character array value
94553acd3b1SBarry Smith 
94653acd3b1SBarry Smith    Output Parameter:
94753acd3b1SBarry Smith +  value - the value to return
94853acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
94953acd3b1SBarry Smith 
95053acd3b1SBarry Smith    Level: intermediate
95153acd3b1SBarry Smith 
95253acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
95353acd3b1SBarry Smith 
95453acd3b1SBarry Smith    See PetscOptionsEList() for when the choices are given in a string array
95553acd3b1SBarry Smith 
95653acd3b1SBarry Smith    To get a listing of all currently specified options,
95788c29154SBarry Smith     see PetscOptionsView() or PetscOptionsGetAll()
95853acd3b1SBarry Smith 
959eabe10d7SBarry Smith    Developer Note: This cannot check for invalid selection because of things like MATAIJ that are not included in the list
960eabe10d7SBarry Smith 
96153acd3b1SBarry Smith    Concepts: options database^list
96253acd3b1SBarry Smith 
96353acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
964acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
96553acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
96653acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
967acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
968a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsEnum()
96953acd3b1SBarry Smith @*/
970a264d7a6SBarry Smith PetscErrorCode  PetscOptionsFList(const char opt[],const char ltext[],const char man[],PetscFunctionList list,const char defaultv[],char value[],size_t len,PetscBool  *set)
97153acd3b1SBarry Smith {
97253acd3b1SBarry Smith   PetscErrorCode ierr;
9733cc1e11dSBarry Smith   PetscOptions   amsopt;
97453acd3b1SBarry Smith 
97553acd3b1SBarry Smith   PetscFunctionBegin;
976b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
977a264d7a6SBarry Smith     ierr = PetscOptionsCreate_Private(opt,ltext,man,OPTION_FLIST,&amsopt);CHKERRQ(ierr);
97864facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
9795b02f95dSBarry Smith     ierr = PetscStrdup(defaultv ? defaultv : "",(char**)&amsopt->data);CHKERRQ(ierr);
9803cc1e11dSBarry Smith     amsopt->flist = list;
9813cc1e11dSBarry Smith   }
98253acd3b1SBarry Smith   ierr = PetscOptionsGetString(PetscOptionsObject.prefix,opt,value,len,set);CHKERRQ(ierr);
98361b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
984140e18c1SBarry Smith     ierr = PetscFunctionListPrintTypes(PetscOptionsObject.comm,stdout,PetscOptionsObject.prefix,opt,ltext,man,list,defaultv);CHKERRQ(ierr);CHKERRQ(ierr);
98553acd3b1SBarry Smith   }
98653acd3b1SBarry Smith   PetscFunctionReturn(0);
98753acd3b1SBarry Smith }
98853acd3b1SBarry Smith 
98953acd3b1SBarry Smith #undef __FUNCT__
99053acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEList"
99153acd3b1SBarry Smith /*@C
99253acd3b1SBarry Smith      PetscOptionsEList - Puts a list of option values that a single one may be selected from
99353acd3b1SBarry Smith 
9943f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
99553acd3b1SBarry Smith 
99653acd3b1SBarry Smith    Input Parameters:
99753acd3b1SBarry Smith +  opt - option name
99853acd3b1SBarry Smith .  ltext - short string that describes the option
99953acd3b1SBarry Smith .  man - manual page with additional information on option
1000a264d7a6SBarry Smith .  list - the possible choices (one of these must be selected, anything else is invalid)
100153acd3b1SBarry Smith .  ntext - number of choices
100253acd3b1SBarry Smith -  defaultv - the default (current) value
100353acd3b1SBarry Smith 
100453acd3b1SBarry Smith    Output Parameter:
100553acd3b1SBarry Smith +  value - the index of the value to return
100653acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
100753acd3b1SBarry Smith 
100853acd3b1SBarry Smith    Level: intermediate
100953acd3b1SBarry Smith 
101053acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
101153acd3b1SBarry Smith 
1012a264d7a6SBarry Smith    See PetscOptionsFList() for when the choices are given in a PetscFunctionList()
101353acd3b1SBarry Smith 
101453acd3b1SBarry Smith    Concepts: options database^list
101553acd3b1SBarry Smith 
101653acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1017acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
101853acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
101953acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1020acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1021a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEnum()
102253acd3b1SBarry Smith @*/
10237087cfbeSBarry 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)
102453acd3b1SBarry Smith {
102553acd3b1SBarry Smith   PetscErrorCode ierr;
102653acd3b1SBarry Smith   PetscInt       i;
10271ae3d29cSBarry Smith   PetscOptions   amsopt;
102853acd3b1SBarry Smith 
102953acd3b1SBarry Smith   PetscFunctionBegin;
10301ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
1031d5649816SBarry Smith     ierr = PetscOptionsCreate_Private(opt,ltext,man,OPTION_ELIST,&amsopt);CHKERRQ(ierr);
103264facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
10335b02f95dSBarry Smith     ierr = PetscStrdup(defaultv ? defaultv : "",(char**)&amsopt->data);CHKERRQ(ierr);
10341ae3d29cSBarry Smith     amsopt->list  = list;
10351ae3d29cSBarry Smith     amsopt->nlist = ntext;
10361ae3d29cSBarry Smith   }
103753acd3b1SBarry Smith   ierr = PetscOptionsGetEList(PetscOptionsObject.prefix,opt,list,ntext,value,set);CHKERRQ(ierr);
103861b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
103953acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%s> (choose one of)",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,defaultv);CHKERRQ(ierr);
104053acd3b1SBarry Smith     for (i=0; i<ntext; i++) {
104153acd3b1SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm," %s",list[i]);CHKERRQ(ierr);
104253acd3b1SBarry Smith     }
10432aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm," (%s)\n",ManSection(man));CHKERRQ(ierr);
104453acd3b1SBarry Smith   }
104553acd3b1SBarry Smith   PetscFunctionReturn(0);
104653acd3b1SBarry Smith }
104753acd3b1SBarry Smith 
104853acd3b1SBarry Smith #undef __FUNCT__
1049acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolGroupBegin"
105053acd3b1SBarry Smith /*@C
1051acfcf0e5SJed Brown      PetscOptionsBoolGroupBegin - First in a series of logical queries on the options database for
1052d5649816SBarry Smith        which at most a single value can be true.
105353acd3b1SBarry Smith 
10543f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
105553acd3b1SBarry Smith 
105653acd3b1SBarry Smith    Input Parameters:
105753acd3b1SBarry Smith +  opt - option name
105853acd3b1SBarry Smith .  text - short string that describes the option
105953acd3b1SBarry Smith -  man - manual page with additional information on option
106053acd3b1SBarry Smith 
106153acd3b1SBarry Smith    Output Parameter:
106253acd3b1SBarry Smith .  flg - whether that option was set or not
106353acd3b1SBarry Smith 
106453acd3b1SBarry Smith    Level: intermediate
106553acd3b1SBarry Smith 
106653acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
106753acd3b1SBarry Smith 
1068acfcf0e5SJed Brown    Must be followed by 0 or more PetscOptionsBoolGroup()s and PetscOptionsBoolGroupEnd()
106953acd3b1SBarry Smith 
107053acd3b1SBarry Smith     Concepts: options database^logical group
107153acd3b1SBarry Smith 
107253acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1073acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
107453acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
107553acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1076acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1077a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
107853acd3b1SBarry Smith @*/
10797087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolGroupBegin(const char opt[],const char text[],const char man[],PetscBool  *flg)
108053acd3b1SBarry Smith {
108153acd3b1SBarry Smith   PetscErrorCode ierr;
10821ae3d29cSBarry Smith   PetscOptions   amsopt;
108353acd3b1SBarry Smith 
108453acd3b1SBarry Smith   PetscFunctionBegin;
10851ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
10867781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1087ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1088a297a907SKarl Rupp 
1089ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
10901ae3d29cSBarry Smith   }
109168b16fdaSBarry Smith   *flg = PETSC_FALSE;
10920298fd71SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,NULL);CHKERRQ(ierr);
109361b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
109453acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  Pick at most one of -------------\n");CHKERRQ(ierr);
10952aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"    -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
109653acd3b1SBarry Smith   }
109753acd3b1SBarry Smith   PetscFunctionReturn(0);
109853acd3b1SBarry Smith }
109953acd3b1SBarry Smith 
110053acd3b1SBarry Smith #undef __FUNCT__
1101acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolGroup"
110253acd3b1SBarry Smith /*@C
1103acfcf0e5SJed Brown      PetscOptionsBoolGroup - One in a series of logical queries on the options database for
1104d5649816SBarry Smith        which at most a single value can be true.
110553acd3b1SBarry Smith 
11063f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
110753acd3b1SBarry Smith 
110853acd3b1SBarry Smith    Input Parameters:
110953acd3b1SBarry Smith +  opt - option name
111053acd3b1SBarry Smith .  text - short string that describes the option
111153acd3b1SBarry Smith -  man - manual page with additional information on option
111253acd3b1SBarry Smith 
111353acd3b1SBarry Smith    Output Parameter:
111453acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
111553acd3b1SBarry Smith 
111653acd3b1SBarry Smith    Level: intermediate
111753acd3b1SBarry Smith 
111853acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
111953acd3b1SBarry Smith 
1120acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin() and preceded a PetscOptionsBoolGroupEnd()
112153acd3b1SBarry Smith 
112253acd3b1SBarry Smith     Concepts: options database^logical group
112353acd3b1SBarry Smith 
112453acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1125acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
112653acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
112753acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1128acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1129a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
113053acd3b1SBarry Smith @*/
11317087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolGroup(const char opt[],const char text[],const char man[],PetscBool  *flg)
113253acd3b1SBarry Smith {
113353acd3b1SBarry Smith   PetscErrorCode ierr;
11341ae3d29cSBarry Smith   PetscOptions   amsopt;
113553acd3b1SBarry Smith 
113653acd3b1SBarry Smith   PetscFunctionBegin;
11371ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
11387781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1139ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1140a297a907SKarl Rupp 
1141ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
11421ae3d29cSBarry Smith   }
114317326d04SJed Brown   *flg = PETSC_FALSE;
11440298fd71SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,NULL);CHKERRQ(ierr);
114561b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
11462aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"    -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
114753acd3b1SBarry Smith   }
114853acd3b1SBarry Smith   PetscFunctionReturn(0);
114953acd3b1SBarry Smith }
115053acd3b1SBarry Smith 
115153acd3b1SBarry Smith #undef __FUNCT__
1152acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolGroupEnd"
115353acd3b1SBarry Smith /*@C
1154acfcf0e5SJed Brown      PetscOptionsBoolGroupEnd - Last in a series of logical queries on the options database for
1155d5649816SBarry Smith        which at most a single value can be true.
115653acd3b1SBarry Smith 
11573f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
115853acd3b1SBarry Smith 
115953acd3b1SBarry Smith    Input Parameters:
116053acd3b1SBarry Smith +  opt - option name
116153acd3b1SBarry Smith .  text - short string that describes the option
116253acd3b1SBarry Smith -  man - manual page with additional information on option
116353acd3b1SBarry Smith 
116453acd3b1SBarry Smith    Output Parameter:
116553acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
116653acd3b1SBarry Smith 
116753acd3b1SBarry Smith    Level: intermediate
116853acd3b1SBarry Smith 
116953acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
117053acd3b1SBarry Smith 
1171acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin()
117253acd3b1SBarry Smith 
117353acd3b1SBarry Smith     Concepts: options database^logical group
117453acd3b1SBarry Smith 
117553acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1176acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
117753acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
117853acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1179acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1180a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
118153acd3b1SBarry Smith @*/
11827087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolGroupEnd(const char opt[],const char text[],const char man[],PetscBool  *flg)
118353acd3b1SBarry Smith {
118453acd3b1SBarry Smith   PetscErrorCode ierr;
11851ae3d29cSBarry Smith   PetscOptions   amsopt;
118653acd3b1SBarry Smith 
118753acd3b1SBarry Smith   PetscFunctionBegin;
11881ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
11897781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1190ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1191a297a907SKarl Rupp 
1192ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
11931ae3d29cSBarry Smith   }
119417326d04SJed Brown   *flg = PETSC_FALSE;
11950298fd71SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,NULL);CHKERRQ(ierr);
119661b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
11972aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"    -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
119853acd3b1SBarry Smith   }
119953acd3b1SBarry Smith   PetscFunctionReturn(0);
120053acd3b1SBarry Smith }
120153acd3b1SBarry Smith 
120253acd3b1SBarry Smith #undef __FUNCT__
1203acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBool"
120453acd3b1SBarry Smith /*@C
1205acfcf0e5SJed Brown    PetscOptionsBool - Determines if a particular option is in the database with a true or false
120653acd3b1SBarry Smith 
12073f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
120853acd3b1SBarry Smith 
120953acd3b1SBarry Smith    Input Parameters:
121053acd3b1SBarry Smith +  opt - option name
121153acd3b1SBarry Smith .  text - short string that describes the option
121253acd3b1SBarry Smith -  man - manual page with additional information on option
121353acd3b1SBarry Smith 
121453acd3b1SBarry Smith    Output Parameter:
121553acd3b1SBarry Smith .  flg - PETSC_TRUE or PETSC_FALSE
121653acd3b1SBarry Smith .  set - PETSC_TRUE if found, else PETSC_FALSE
121753acd3b1SBarry Smith 
121853acd3b1SBarry Smith    Level: beginner
121953acd3b1SBarry Smith 
122053acd3b1SBarry Smith    Concepts: options database^logical
122153acd3b1SBarry Smith 
122253acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
122353acd3b1SBarry Smith 
122453acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
1225acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
1226acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
122753acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
122853acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1229acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1230a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
123153acd3b1SBarry Smith @*/
12327087cfbeSBarry Smith PetscErrorCode  PetscOptionsBool(const char opt[],const char text[],const char man[],PetscBool deflt,PetscBool  *flg,PetscBool  *set)
123353acd3b1SBarry Smith {
123453acd3b1SBarry Smith   PetscErrorCode ierr;
1235ace3abfcSBarry Smith   PetscBool      iset;
1236aee2cecaSBarry Smith   PetscOptions   amsopt;
123753acd3b1SBarry Smith 
123853acd3b1SBarry Smith   PetscFunctionBegin;
1239b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
12407781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1241ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1242a297a907SKarl Rupp 
1243ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = deflt;
1244af6d86caSBarry Smith   }
1245acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,&iset);CHKERRQ(ierr);
124653acd3b1SBarry Smith   if (!iset) {
124753acd3b1SBarry Smith     if (flg) *flg = deflt;
124853acd3b1SBarry Smith   }
124953acd3b1SBarry Smith   if (set) *set = iset;
125061b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
1251ace3abfcSBarry Smith     const char *v = PetscBools[deflt];
12522aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s: <%s> %s (%s)\n",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,v,text,ManSection(man));CHKERRQ(ierr);
125353acd3b1SBarry Smith   }
125453acd3b1SBarry Smith   PetscFunctionReturn(0);
125553acd3b1SBarry Smith }
125653acd3b1SBarry Smith 
125753acd3b1SBarry Smith #undef __FUNCT__
125853acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsRealArray"
125953acd3b1SBarry Smith /*@C
126053acd3b1SBarry Smith    PetscOptionsRealArray - Gets an array of double values for a particular
126153acd3b1SBarry Smith    option in the database. The values must be separated with commas with
126253acd3b1SBarry Smith    no intervening spaces.
126353acd3b1SBarry Smith 
12643f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
126553acd3b1SBarry Smith 
126653acd3b1SBarry Smith    Input Parameters:
126753acd3b1SBarry Smith +  opt - the option one is seeking
126853acd3b1SBarry Smith .  text - short string describing option
126953acd3b1SBarry Smith .  man - manual page for option
127053acd3b1SBarry Smith -  nmax - maximum number of values
127153acd3b1SBarry Smith 
127253acd3b1SBarry Smith    Output Parameter:
127353acd3b1SBarry Smith +  value - location to copy values
127453acd3b1SBarry Smith .  nmax - actual number of values found
127553acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
127653acd3b1SBarry Smith 
127753acd3b1SBarry Smith    Level: beginner
127853acd3b1SBarry Smith 
127953acd3b1SBarry Smith    Notes:
128053acd3b1SBarry Smith    The user should pass in an array of doubles
128153acd3b1SBarry Smith 
128253acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
128353acd3b1SBarry Smith 
128453acd3b1SBarry Smith    Concepts: options database^array of strings
128553acd3b1SBarry Smith 
128653acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1287acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
128853acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
128953acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1290acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1291a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
129253acd3b1SBarry Smith @*/
12937087cfbeSBarry Smith PetscErrorCode  PetscOptionsRealArray(const char opt[],const char text[],const char man[],PetscReal value[],PetscInt *n,PetscBool  *set)
129453acd3b1SBarry Smith {
129553acd3b1SBarry Smith   PetscErrorCode ierr;
129653acd3b1SBarry Smith   PetscInt       i;
1297e26ddf31SBarry Smith   PetscOptions   amsopt;
129853acd3b1SBarry Smith 
129953acd3b1SBarry Smith   PetscFunctionBegin;
1300b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
1301e26ddf31SBarry Smith     PetscReal *vals;
1302e26ddf31SBarry Smith 
1303e26ddf31SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_REAL_ARRAY,&amsopt);CHKERRQ(ierr);
13046e02237eSPeter Brune     ierr = PetscMalloc1((*n),(PetscReal**)&amsopt->data);CHKERRQ(ierr);
1305e26ddf31SBarry Smith     vals = (PetscReal*)amsopt->data;
1306e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1307e26ddf31SBarry Smith     amsopt->arraylength = *n;
1308e26ddf31SBarry Smith   }
130953acd3b1SBarry Smith   ierr = PetscOptionsGetRealArray(PetscOptionsObject.prefix,opt,value,n,set);CHKERRQ(ierr);
131061b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
1311a83599f4SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%G",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,value[0]);CHKERRQ(ierr);
131253acd3b1SBarry Smith     for (i=1; i<*n; i++) {
1313a83599f4SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,",%G",value[i]);CHKERRQ(ierr);
131453acd3b1SBarry Smith     }
13152aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
131653acd3b1SBarry Smith   }
131753acd3b1SBarry Smith   PetscFunctionReturn(0);
131853acd3b1SBarry Smith }
131953acd3b1SBarry Smith 
132053acd3b1SBarry Smith 
132153acd3b1SBarry Smith #undef __FUNCT__
132253acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsIntArray"
132353acd3b1SBarry Smith /*@C
132453acd3b1SBarry Smith    PetscOptionsIntArray - Gets an array of integers for a particular
1325b32a342fSShri Abhyankar    option in the database.
132653acd3b1SBarry Smith 
13273f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
132853acd3b1SBarry Smith 
132953acd3b1SBarry Smith    Input Parameters:
133053acd3b1SBarry Smith +  opt - the option one is seeking
133153acd3b1SBarry Smith .  text - short string describing option
133253acd3b1SBarry Smith .  man - manual page for option
1333f8a50e2bSBarry Smith -  n - maximum number of values
133453acd3b1SBarry Smith 
133553acd3b1SBarry Smith    Output Parameter:
133653acd3b1SBarry Smith +  value - location to copy values
1337f8a50e2bSBarry Smith .  n - actual number of values found
133853acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
133953acd3b1SBarry Smith 
134053acd3b1SBarry Smith    Level: beginner
134153acd3b1SBarry Smith 
134253acd3b1SBarry Smith    Notes:
1343b32a342fSShri Abhyankar    The array can be passed as
1344b32a342fSShri Abhyankar    a comma seperated list:                                 0,1,2,3,4,5,6,7
13450fd488f5SShri Abhyankar    a range (start-end+1):                                  0-8
13460fd488f5SShri Abhyankar    a range with given increment (start-end+1:inc):         0-7:2
13470fd488f5SShri Abhyankar    a combination of values and ranges seperated by commas: 0,1-8,8-15:2
1348b32a342fSShri Abhyankar 
1349b32a342fSShri Abhyankar    There must be no intervening spaces between the values.
135053acd3b1SBarry Smith 
135153acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
135253acd3b1SBarry Smith 
1353b32a342fSShri Abhyankar    Concepts: options database^array of ints
135453acd3b1SBarry Smith 
135553acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1356acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
135753acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
135853acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1359acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1360a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsRealArray()
136153acd3b1SBarry Smith @*/
13627087cfbeSBarry Smith PetscErrorCode  PetscOptionsIntArray(const char opt[],const char text[],const char man[],PetscInt value[],PetscInt *n,PetscBool  *set)
136353acd3b1SBarry Smith {
136453acd3b1SBarry Smith   PetscErrorCode ierr;
136553acd3b1SBarry Smith   PetscInt       i;
1366e26ddf31SBarry Smith   PetscOptions   amsopt;
136753acd3b1SBarry Smith 
136853acd3b1SBarry Smith   PetscFunctionBegin;
1369b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
1370e26ddf31SBarry Smith     PetscInt *vals;
1371e26ddf31SBarry Smith 
1372e26ddf31SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_INT_ARRAY,&amsopt);CHKERRQ(ierr);
13736e02237eSPeter Brune     ierr = PetscMalloc1((*n),(PetscInt**)&amsopt->data);CHKERRQ(ierr);
1374e26ddf31SBarry Smith     vals = (PetscInt*)amsopt->data;
1375e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1376e26ddf31SBarry Smith     amsopt->arraylength = *n;
1377e26ddf31SBarry Smith   }
137853acd3b1SBarry Smith   ierr = PetscOptionsGetIntArray(PetscOptionsObject.prefix,opt,value,n,set);CHKERRQ(ierr);
137961b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
138053acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%d",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,value[0]);CHKERRQ(ierr);
138153acd3b1SBarry Smith     for (i=1; i<*n; i++) {
138253acd3b1SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,",%d",value[i]);CHKERRQ(ierr);
138353acd3b1SBarry Smith     }
13842aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
138553acd3b1SBarry Smith   }
138653acd3b1SBarry Smith   PetscFunctionReturn(0);
138753acd3b1SBarry Smith }
138853acd3b1SBarry Smith 
138953acd3b1SBarry Smith #undef __FUNCT__
139053acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsStringArray"
139153acd3b1SBarry Smith /*@C
139253acd3b1SBarry Smith    PetscOptionsStringArray - Gets an array of string values for a particular
139353acd3b1SBarry Smith    option in the database. The values must be separated with commas with
139453acd3b1SBarry Smith    no intervening spaces.
139553acd3b1SBarry Smith 
13963f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
139753acd3b1SBarry Smith 
139853acd3b1SBarry Smith    Input Parameters:
139953acd3b1SBarry Smith +  opt - the option one is seeking
140053acd3b1SBarry Smith .  text - short string describing option
140153acd3b1SBarry Smith .  man - manual page for option
140253acd3b1SBarry Smith -  nmax - maximum number of strings
140353acd3b1SBarry Smith 
140453acd3b1SBarry Smith    Output Parameter:
140553acd3b1SBarry Smith +  value - location to copy strings
140653acd3b1SBarry Smith .  nmax - actual number of strings found
140753acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
140853acd3b1SBarry Smith 
140953acd3b1SBarry Smith    Level: beginner
141053acd3b1SBarry Smith 
141153acd3b1SBarry Smith    Notes:
141253acd3b1SBarry Smith    The user should pass in an array of pointers to char, to hold all the
141353acd3b1SBarry Smith    strings returned by this function.
141453acd3b1SBarry Smith 
141553acd3b1SBarry Smith    The user is responsible for deallocating the strings that are
141653acd3b1SBarry Smith    returned. The Fortran interface for this routine is not supported.
141753acd3b1SBarry Smith 
141853acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
141953acd3b1SBarry Smith 
142053acd3b1SBarry Smith    Concepts: options database^array of strings
142153acd3b1SBarry Smith 
142253acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1423acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
142453acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
142553acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1426acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1427a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
142853acd3b1SBarry Smith @*/
14297087cfbeSBarry Smith PetscErrorCode  PetscOptionsStringArray(const char opt[],const char text[],const char man[],char *value[],PetscInt *nmax,PetscBool  *set)
143053acd3b1SBarry Smith {
143153acd3b1SBarry Smith   PetscErrorCode ierr;
14321ae3d29cSBarry Smith   PetscOptions   amsopt;
143353acd3b1SBarry Smith 
143453acd3b1SBarry Smith   PetscFunctionBegin;
14351ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
14361ae3d29cSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_STRING_ARRAY,&amsopt);CHKERRQ(ierr);
14376e02237eSPeter Brune     ierr = PetscMalloc1((*nmax),(char**)&amsopt->data);CHKERRQ(ierr);
1438a297a907SKarl Rupp 
14391ae3d29cSBarry Smith     amsopt->arraylength = *nmax;
14401ae3d29cSBarry Smith   }
144153acd3b1SBarry Smith   ierr = PetscOptionsGetStringArray(PetscOptionsObject.prefix,opt,value,nmax,set);CHKERRQ(ierr);
144261b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
14432aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <string1,string2,...>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
144453acd3b1SBarry Smith   }
144553acd3b1SBarry Smith   PetscFunctionReturn(0);
144653acd3b1SBarry Smith }
144753acd3b1SBarry Smith 
1448e2446a98SMatthew Knepley #undef __FUNCT__
1449acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolArray"
1450e2446a98SMatthew Knepley /*@C
1451acfcf0e5SJed Brown    PetscOptionsBoolArray - Gets an array of logical values (true or false) for a particular
1452e2446a98SMatthew Knepley    option in the database. The values must be separated with commas with
1453e2446a98SMatthew Knepley    no intervening spaces.
1454e2446a98SMatthew Knepley 
14553f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
1456e2446a98SMatthew Knepley 
1457e2446a98SMatthew Knepley    Input Parameters:
1458e2446a98SMatthew Knepley +  opt - the option one is seeking
1459e2446a98SMatthew Knepley .  text - short string describing option
1460e2446a98SMatthew Knepley .  man - manual page for option
1461e2446a98SMatthew Knepley -  nmax - maximum number of values
1462e2446a98SMatthew Knepley 
1463e2446a98SMatthew Knepley    Output Parameter:
1464e2446a98SMatthew Knepley +  value - location to copy values
1465e2446a98SMatthew Knepley .  nmax - actual number of values found
1466e2446a98SMatthew Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
1467e2446a98SMatthew Knepley 
1468e2446a98SMatthew Knepley    Level: beginner
1469e2446a98SMatthew Knepley 
1470e2446a98SMatthew Knepley    Notes:
1471e2446a98SMatthew Knepley    The user should pass in an array of doubles
1472e2446a98SMatthew Knepley 
1473e2446a98SMatthew Knepley    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
1474e2446a98SMatthew Knepley 
1475e2446a98SMatthew Knepley    Concepts: options database^array of strings
1476e2446a98SMatthew Knepley 
1477e2446a98SMatthew Knepley .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1478acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
1479e2446a98SMatthew Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1480e2446a98SMatthew Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1481acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1482a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
1483e2446a98SMatthew Knepley @*/
14847087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolArray(const char opt[],const char text[],const char man[],PetscBool value[],PetscInt *n,PetscBool *set)
1485e2446a98SMatthew Knepley {
1486e2446a98SMatthew Knepley   PetscErrorCode ierr;
1487e2446a98SMatthew Knepley   PetscInt       i;
14881ae3d29cSBarry Smith   PetscOptions   amsopt;
1489e2446a98SMatthew Knepley 
1490e2446a98SMatthew Knepley   PetscFunctionBegin;
14911ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
1492ace3abfcSBarry Smith     PetscBool *vals;
14931ae3d29cSBarry Smith 
14947781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL_ARRAY,&amsopt);CHKERRQ(ierr);
14956e02237eSPeter Brune     ierr = PetscMalloc1((*n),(PetscBool**)&amsopt->data);CHKERRQ(ierr);
1496ace3abfcSBarry Smith     vals = (PetscBool*)amsopt->data;
14971ae3d29cSBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
14981ae3d29cSBarry Smith     amsopt->arraylength = *n;
14991ae3d29cSBarry Smith   }
1500acfcf0e5SJed Brown   ierr = PetscOptionsGetBoolArray(PetscOptionsObject.prefix,opt,value,n,set);CHKERRQ(ierr);
1501e2446a98SMatthew Knepley   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
1502e2446a98SMatthew Knepley     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%d",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,value[0]);CHKERRQ(ierr);
1503e2446a98SMatthew Knepley     for (i=1; i<*n; i++) {
1504e2446a98SMatthew Knepley       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,",%d",value[i]);CHKERRQ(ierr);
1505e2446a98SMatthew Knepley     }
15062aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
1507e2446a98SMatthew Knepley   }
1508e2446a98SMatthew Knepley   PetscFunctionReturn(0);
1509e2446a98SMatthew Knepley }
1510e2446a98SMatthew Knepley 
15118cc676e6SMatthew G Knepley #undef __FUNCT__
15128cc676e6SMatthew G Knepley #define __FUNCT__ "PetscOptionsViewer"
15138cc676e6SMatthew G Knepley /*@C
15148cc676e6SMatthew G Knepley    PetscOptionsInt - Gets a viewer appropriate for the type indicated by the user
15158cc676e6SMatthew G Knepley 
15168cc676e6SMatthew G Knepley    Logically Collective on the communicator passed in PetscOptionsBegin()
15178cc676e6SMatthew G Knepley 
15188cc676e6SMatthew G Knepley    Input Parameters:
15198cc676e6SMatthew G Knepley +  opt - option name
15208cc676e6SMatthew G Knepley .  text - short string that describes the option
15218cc676e6SMatthew G Knepley -  man - manual page with additional information on option
15228cc676e6SMatthew G Knepley 
15238cc676e6SMatthew G Knepley    Output Parameter:
15248cc676e6SMatthew G Knepley +  viewer - the viewer
15258cc676e6SMatthew G Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
15268cc676e6SMatthew G Knepley 
15278cc676e6SMatthew G Knepley    Level: beginner
15288cc676e6SMatthew G Knepley 
15298cc676e6SMatthew G Knepley    Concepts: options database^has int
15308cc676e6SMatthew G Knepley 
15318cc676e6SMatthew G Knepley    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
15328cc676e6SMatthew G Knepley If no value is provided ascii:stdout is used
15338cc676e6SMatthew 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
15348cc676e6SMatthew G Knepley $                                     about the object to standard out
15358cc676e6SMatthew G Knepley $       binary[:filename]   defaults to binaryoutput
15368cc676e6SMatthew G Knepley $       draw
15378cc676e6SMatthew G Knepley $       socket[:port]    defaults to the standard output port
15388cc676e6SMatthew G Knepley 
1539cffb1e40SBarry Smith    Use PetscRestoreViewerDestroy() after using the viewer, otherwise a memory leak will occur
15408cc676e6SMatthew G Knepley 
15418cc676e6SMatthew G Knepley .seealso: PetscOptionsGetViewer(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
15428cc676e6SMatthew G Knepley           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
15438cc676e6SMatthew G Knepley           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
15448cc676e6SMatthew G Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
15458cc676e6SMatthew G Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
15468cc676e6SMatthew G Knepley           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1547a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
15488cc676e6SMatthew G Knepley @*/
1549cffb1e40SBarry Smith PetscErrorCode  PetscOptionsViewer(const char opt[],const char text[],const char man[],PetscViewer *viewer,PetscViewerFormat *format,PetscBool  *set)
15508cc676e6SMatthew G Knepley {
15518cc676e6SMatthew G Knepley   PetscErrorCode ierr;
15528cc676e6SMatthew G Knepley   PetscOptions   amsopt;
15538cc676e6SMatthew G Knepley 
15548cc676e6SMatthew G Knepley   PetscFunctionBegin;
15558cc676e6SMatthew G Knepley   if (!PetscOptionsPublishCount) {
15568cc676e6SMatthew G Knepley     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
155764facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
15585b02f95dSBarry Smith     ierr = PetscStrdup("",(char**)&amsopt->data);CHKERRQ(ierr);
15598cc676e6SMatthew G Knepley   }
1560cffb1e40SBarry Smith   ierr = PetscOptionsGetViewer(PetscOptionsObject.comm,PetscOptionsObject.prefix,opt,viewer,format,set);CHKERRQ(ierr);
15618cc676e6SMatthew G Knepley   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
15628cc676e6SMatthew G Knepley     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%s>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,"",text,ManSection(man));CHKERRQ(ierr);
15638cc676e6SMatthew G Knepley   }
15648cc676e6SMatthew G Knepley   PetscFunctionReturn(0);
15658cc676e6SMatthew G Knepley }
15668cc676e6SMatthew G Knepley 
156753acd3b1SBarry Smith 
156853acd3b1SBarry Smith #undef __FUNCT__
156953acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsHead"
157053acd3b1SBarry Smith /*@C
1571b52f573bSBarry Smith      PetscOptionsHead - Puts a heading before listing any more published options. Used, for example,
157253acd3b1SBarry Smith             in KSPSetFromOptions_GMRES().
157353acd3b1SBarry Smith 
15743f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
157553acd3b1SBarry Smith 
157653acd3b1SBarry Smith    Input Parameter:
157753acd3b1SBarry Smith .   head - the heading text
157853acd3b1SBarry Smith 
157953acd3b1SBarry Smith 
158053acd3b1SBarry Smith    Level: intermediate
158153acd3b1SBarry Smith 
158253acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
158353acd3b1SBarry Smith 
1584b52f573bSBarry Smith           Can be followed by a call to PetscOptionsTail() in the same function.
158553acd3b1SBarry Smith 
158653acd3b1SBarry Smith    Concepts: options database^subheading
158753acd3b1SBarry Smith 
158853acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1589acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
159053acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
159153acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1592acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1593a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
159453acd3b1SBarry Smith @*/
15957087cfbeSBarry Smith PetscErrorCode  PetscOptionsHead(const char head[])
159653acd3b1SBarry Smith {
159753acd3b1SBarry Smith   PetscErrorCode ierr;
159853acd3b1SBarry Smith 
159953acd3b1SBarry Smith   PetscFunctionBegin;
160061b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
160153acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  %s\n",head);CHKERRQ(ierr);
160253acd3b1SBarry Smith   }
160353acd3b1SBarry Smith   PetscFunctionReturn(0);
160453acd3b1SBarry Smith }
160553acd3b1SBarry Smith 
160653acd3b1SBarry Smith 
160753acd3b1SBarry Smith 
160853acd3b1SBarry Smith 
160953acd3b1SBarry Smith 
161053acd3b1SBarry Smith 
1611