xref: /petsc/src/sys/objects/aoptions.c (revision aca04275f8154a79110c83ebd0913ae08fc5b083)
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 
149*aca04275SShao-Ching Huang #undef __FUNCT__
150*aca04275SShao-Ching Huang #define __FUNCT__ "PetscStrdup"
1515b02f95dSBarry Smith /*
1525b02f95dSBarry Smith     This is needed because certain strings may be freed by SAWs, hence we cannot use PetscStrallocpy()
1535b02f95dSBarry Smith */
1545b02f95dSBarry Smith static PetscErrorCode  PetscStrdup(const char s[],char *t[])
1555b02f95dSBarry Smith {
1565b02f95dSBarry Smith   PetscErrorCode ierr;
1575b02f95dSBarry Smith   size_t         len;
1585b02f95dSBarry Smith   char           *tmp = 0;
1595b02f95dSBarry Smith 
1605b02f95dSBarry Smith   PetscFunctionBegin;
1615b02f95dSBarry Smith   if (s) {
1625b02f95dSBarry Smith     ierr = PetscStrlen(s,&len);CHKERRQ(ierr);
163ee48884fSBarry Smith     tmp = (char*) malloc((len+1)*sizeof(char*));
1645b02f95dSBarry Smith     if (!tmp) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"No memory to duplicate string");
1655b02f95dSBarry Smith     ierr = PetscStrcpy(tmp,s);CHKERRQ(ierr);
1665b02f95dSBarry Smith   }
1675b02f95dSBarry Smith   *t = tmp;
1685b02f95dSBarry Smith   PetscFunctionReturn(0);
1695b02f95dSBarry Smith }
1705b02f95dSBarry Smith 
1715b02f95dSBarry Smith 
172aee2cecaSBarry Smith #undef __FUNCT__
173aee2cecaSBarry Smith #define __FUNCT__ "PetscOptionsGetFromTextInput"
174aee2cecaSBarry Smith /*
1753cc1e11dSBarry Smith     PetscOptionsGetFromTextInput - Presents all the PETSc Options processed by the program so the user may change them at runtime
176aee2cecaSBarry Smith 
177aee2cecaSBarry Smith     Notes: this isn't really practical, it is just to demonstrate the principle
178aee2cecaSBarry Smith 
1797781c08eSBarry Smith     A carriage return indicates no change from the default; but this like -ksp_monitor <stdout>  the default is actually not stdout the default
1807781c08eSBarry Smith     is to do nothing so to get it to use stdout you need to type stdout. This is kind of bug?
1817781c08eSBarry Smith 
182aee2cecaSBarry Smith     Bugs:
1837781c08eSBarry Smith +    All processes must traverse through the exact same set of option queries due to the call to PetscScanString()
1843cc1e11dSBarry Smith .    Internal strings have arbitrary length and string copies are not checked that they fit into string space
185aee2cecaSBarry Smith -    Only works for PetscInt == int, PetscReal == double etc
186aee2cecaSBarry Smith 
1873cc1e11dSBarry Smith     Developer Notes: Normally the GUI that presents the options the user and retrieves the values would be running in a different
1883cc1e11dSBarry Smith      address space and communicating with the PETSc program
1893cc1e11dSBarry Smith 
190aee2cecaSBarry Smith */
191aee2cecaSBarry Smith PetscErrorCode PetscOptionsGetFromTextInput()
1926356e834SBarry Smith {
1936356e834SBarry Smith   PetscErrorCode ierr;
1946356e834SBarry Smith   PetscOptions   next = PetscOptionsObject.next;
1956356e834SBarry Smith   char           str[512];
1967781c08eSBarry Smith   PetscBool      bid;
197a4404d99SBarry Smith   PetscReal      ir,*valr;
198330cf3c9SBarry Smith   PetscInt       *vald;
199330cf3c9SBarry Smith   size_t         i;
2006356e834SBarry Smith 
201e26ddf31SBarry Smith   ierr = (*PetscPrintf)(PETSC_COMM_WORLD,"%s -------------------------------------------------\n",PetscOptionsObject.title);CHKERRQ(ierr);
2026356e834SBarry Smith   while (next) {
2036356e834SBarry Smith     switch (next->type) {
2046356e834SBarry Smith     case OPTION_HEAD:
2056356e834SBarry Smith       break;
206e26ddf31SBarry Smith     case OPTION_INT_ARRAY:
207e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",next->option+1);CHKERRQ(ierr);
208e26ddf31SBarry Smith       vald = (PetscInt*) next->data;
209e26ddf31SBarry Smith       for (i=0; i<next->arraylength; i++) {
210e26ddf31SBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"%d",vald[i]);CHKERRQ(ierr);
211e26ddf31SBarry Smith         if (i < next->arraylength-1) {
212e26ddf31SBarry Smith           ierr = PetscPrintf(PETSC_COMM_WORLD,",");CHKERRQ(ierr);
213e26ddf31SBarry Smith         }
214e26ddf31SBarry Smith       }
215e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,">: %s (%s) ",next->text,next->man);CHKERRQ(ierr);
216e26ddf31SBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
217e26ddf31SBarry Smith       if (str[0]) {
218e26ddf31SBarry Smith         PetscToken token;
219e26ddf31SBarry Smith         PetscInt   n=0,nmax = next->arraylength,*dvalue = (PetscInt*)next->data,start,end;
220e26ddf31SBarry Smith         size_t     len;
221e26ddf31SBarry Smith         char       *value;
222ace3abfcSBarry Smith         PetscBool  foundrange;
223e26ddf31SBarry Smith 
224e26ddf31SBarry Smith         next->set = PETSC_TRUE;
225e26ddf31SBarry Smith         value     = str;
226e26ddf31SBarry Smith         ierr      = PetscTokenCreate(value,',',&token);CHKERRQ(ierr);
227e26ddf31SBarry Smith         ierr      = PetscTokenFind(token,&value);CHKERRQ(ierr);
228e26ddf31SBarry Smith         while (n < nmax) {
229e26ddf31SBarry Smith           if (!value) break;
230e26ddf31SBarry Smith 
231e26ddf31SBarry Smith           /* look for form  d-D where d and D are integers */
232e26ddf31SBarry Smith           foundrange = PETSC_FALSE;
233e26ddf31SBarry Smith           ierr       = PetscStrlen(value,&len);CHKERRQ(ierr);
234e26ddf31SBarry Smith           if (value[0] == '-') i=2;
235e26ddf31SBarry Smith           else i=1;
236330cf3c9SBarry Smith           for (;i<len; i++) {
237e26ddf31SBarry Smith             if (value[i] == '-') {
238e32f2f54SBarry Smith               if (i == len-1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Error in %D-th array entry %s\n",n,value);
239e26ddf31SBarry Smith               value[i] = 0;
240cfbddea1SSatish Balay               ierr     = PetscOptionsStringToInt(value,&start);CHKERRQ(ierr);
241cfbddea1SSatish Balay               ierr     = PetscOptionsStringToInt(value+i+1,&end);CHKERRQ(ierr);
242e32f2f54SBarry Smith               if (end <= start) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_USER,"Error in %D-th array entry, %s-%s cannot have decreasing list",n,value,value+i+1);
243e32f2f54SBarry Smith               if (n + end - start - 1 >= nmax) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_USER,"Error in %D-th array entry, not enough space in left in array (%D) to contain entire range from %D to %D",n,nmax-n,start,end);
244e26ddf31SBarry Smith               for (; start<end; start++) {
245e26ddf31SBarry Smith                 *dvalue = start; dvalue++;n++;
246e26ddf31SBarry Smith               }
247e26ddf31SBarry Smith               foundrange = PETSC_TRUE;
248e26ddf31SBarry Smith               break;
249e26ddf31SBarry Smith             }
250e26ddf31SBarry Smith           }
251e26ddf31SBarry Smith           if (!foundrange) {
252cfbddea1SSatish Balay             ierr = PetscOptionsStringToInt(value,dvalue);CHKERRQ(ierr);
253e26ddf31SBarry Smith             dvalue++;
254e26ddf31SBarry Smith             n++;
255e26ddf31SBarry Smith           }
256e26ddf31SBarry Smith           ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
257e26ddf31SBarry Smith         }
2588c74ee41SBarry Smith         ierr = PetscTokenDestroy(&token);CHKERRQ(ierr);
259e26ddf31SBarry Smith       }
260e26ddf31SBarry Smith       break;
261e26ddf31SBarry Smith     case OPTION_REAL_ARRAY:
262e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",next->option+1);CHKERRQ(ierr);
263e26ddf31SBarry Smith       valr = (PetscReal*) next->data;
264e26ddf31SBarry Smith       for (i=0; i<next->arraylength; i++) {
265e26ddf31SBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"%g",valr[i]);CHKERRQ(ierr);
266e26ddf31SBarry Smith         if (i < next->arraylength-1) {
267e26ddf31SBarry Smith           ierr = PetscPrintf(PETSC_COMM_WORLD,",");CHKERRQ(ierr);
268e26ddf31SBarry Smith         }
269e26ddf31SBarry Smith       }
270e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,">: %s (%s) ",next->text,next->man);CHKERRQ(ierr);
271e26ddf31SBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
272e26ddf31SBarry Smith       if (str[0]) {
273e26ddf31SBarry Smith         PetscToken token;
274e26ddf31SBarry Smith         PetscInt   n = 0,nmax = next->arraylength;
275e26ddf31SBarry Smith         PetscReal  *dvalue = (PetscReal*)next->data;
276e26ddf31SBarry Smith         char       *value;
277e26ddf31SBarry Smith 
278e26ddf31SBarry Smith         next->set = PETSC_TRUE;
279e26ddf31SBarry Smith         value     = str;
280e26ddf31SBarry Smith         ierr      = PetscTokenCreate(value,',',&token);CHKERRQ(ierr);
281e26ddf31SBarry Smith         ierr      = PetscTokenFind(token,&value);CHKERRQ(ierr);
282e26ddf31SBarry Smith         while (n < nmax) {
283e26ddf31SBarry Smith           if (!value) break;
284cfbddea1SSatish Balay           ierr = PetscOptionsStringToReal(value,dvalue);CHKERRQ(ierr);
285e26ddf31SBarry Smith           dvalue++;
286e26ddf31SBarry Smith           n++;
287e26ddf31SBarry Smith           ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
288e26ddf31SBarry Smith         }
2898c74ee41SBarry Smith         ierr = PetscTokenDestroy(&token);CHKERRQ(ierr);
290e26ddf31SBarry Smith       }
291e26ddf31SBarry Smith       break;
2926356e834SBarry Smith     case OPTION_INT:
293e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <%d>: %s (%s) ",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",next->option+1,*(int*)next->data,next->text,next->man);CHKERRQ(ierr);
2943fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
2953fc1eb6aSBarry Smith       if (str[0]) {
296d25d7f95SJed Brown #if defined(PETSC_SIZEOF_LONG_LONG)
297d25d7f95SJed Brown         long long lid;
298d25d7f95SJed Brown         sscanf(str,"%lld",&lid);
299d25d7f95SJed Brown         if (lid > PETSC_MAX_INT || lid < PETSC_MIN_INT) SETERRQ3(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Argument: -%s%s %lld",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",next->option+1,lid);
300c272547aSJed Brown #else
301d25d7f95SJed Brown         long  lid;
302d25d7f95SJed Brown         sscanf(str,"%ld",&lid);
303d25d7f95SJed Brown         if (lid > PETSC_MAX_INT || lid < PETSC_MIN_INT) SETERRQ3(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Argument: -%s%s %ld",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",next->option+1,lid);
304c272547aSJed Brown #endif
305a297a907SKarl Rupp 
306d25d7f95SJed Brown         next->set = PETSC_TRUE;
307d25d7f95SJed Brown         *((PetscInt*)next->data) = (PetscInt)lid;
308aee2cecaSBarry Smith       }
309aee2cecaSBarry Smith       break;
310aee2cecaSBarry Smith     case OPTION_REAL:
311e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <%g>: %s (%s) ",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",next->option+1,*(double*)next->data,next->text,next->man);CHKERRQ(ierr);
3123fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3133fc1eb6aSBarry Smith       if (str[0]) {
314ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
315a4404d99SBarry Smith         sscanf(str,"%e",&ir);
316ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL_DOUBLE)
317aee2cecaSBarry Smith         sscanf(str,"%le",&ir);
318ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL___FLOAT128)
319d9822059SBarry Smith         ir = strtoflt128(str,0);
320d9822059SBarry Smith #else
321513dbe71SLisandro Dalcin         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Unknown scalar type");
322a4404d99SBarry Smith #endif
323aee2cecaSBarry Smith         next->set                 = PETSC_TRUE;
324aee2cecaSBarry Smith         *((PetscReal*)next->data) = ir;
325aee2cecaSBarry Smith       }
326aee2cecaSBarry Smith       break;
3277781c08eSBarry Smith     case OPTION_BOOL:
3287781c08eSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <%s>: %s (%s) ",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",next->option+1,*(PetscBool*)next->data ? "true": "false",next->text,next->man);CHKERRQ(ierr);
3297781c08eSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3307781c08eSBarry Smith       if (str[0]) {
3317781c08eSBarry Smith         ierr = PetscOptionsStringToBool(str,&bid);CHKERRQ(ierr);
3327781c08eSBarry Smith         next->set = PETSC_TRUE;
3337781c08eSBarry Smith         *((PetscBool*)next->data) = bid;
3347781c08eSBarry Smith       }
3357781c08eSBarry Smith       break;
336aee2cecaSBarry Smith     case OPTION_STRING:
337e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <%s>: %s (%s) ",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",next->option+1,(char*)next->data,next->text,next->man);CHKERRQ(ierr);
3383fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3393fc1eb6aSBarry Smith       if (str[0]) {
340aee2cecaSBarry Smith         next->set = PETSC_TRUE;
34164facd6cSBarry Smith         /* must use system malloc since SAWs may free this */
3425b02f95dSBarry Smith         ierr = PetscStrdup(str,(char**)&next->data);CHKERRQ(ierr);
3436356e834SBarry Smith       }
3446356e834SBarry Smith       break;
345a264d7a6SBarry Smith     case OPTION_FLIST:
346140e18c1SBarry Smith       ierr = PetscFunctionListPrintTypes(PETSC_COMM_WORLD,stdout,PetscOptionsObject.prefix,next->option,next->text,next->man,next->flist,(char*)next->data);CHKERRQ(ierr);
3473cc1e11dSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3483cc1e11dSBarry Smith       if (str[0]) {
3493cc1e11dSBarry Smith         PetscOptionsObject.changedmethod = PETSC_TRUE;
3503cc1e11dSBarry Smith         next->set = PETSC_TRUE;
35164facd6cSBarry Smith         /* must use system malloc since SAWs may free this */
3525b02f95dSBarry Smith         ierr = PetscStrdup(str,(char**)&next->data);CHKERRQ(ierr);
3533cc1e11dSBarry Smith       }
3543cc1e11dSBarry Smith       break;
355b432afdaSMatthew Knepley     default:
356b432afdaSMatthew Knepley       break;
3576356e834SBarry Smith     }
3586356e834SBarry Smith     next = next->next;
3596356e834SBarry Smith   }
3606356e834SBarry Smith   PetscFunctionReturn(0);
3616356e834SBarry Smith }
3626356e834SBarry Smith 
363e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
364e04113cfSBarry Smith #include <petscviewersaws.h>
365d5649816SBarry Smith 
366d5649816SBarry Smith static int count = 0;
367d5649816SBarry Smith 
368b3506946SBarry Smith #undef __FUNCT__
369e04113cfSBarry Smith #define __FUNCT__ "PetscOptionsSAWsDestroy"
370e04113cfSBarry Smith PetscErrorCode PetscOptionsSAWsDestroy(void)
371d5649816SBarry Smith {
3722657e9d9SBarry Smith   PetscFunctionBegin;
373d5649816SBarry Smith   PetscFunctionReturn(0);
374d5649816SBarry Smith }
375d5649816SBarry Smith 
376d5649816SBarry Smith #undef __FUNCT__
3777781c08eSBarry Smith #define __FUNCT__ "PetscOptionsSAWsInput"
378b3506946SBarry Smith /*
3797781c08eSBarry Smith     PetscOptionsSAWsInput - Presents all the PETSc Options processed by the program so the user may change them at runtime using the SAWs
380b3506946SBarry Smith 
381b3506946SBarry Smith     Bugs:
382b3506946SBarry Smith +    All processes must traverse through the exact same set of option queries do to the call to PetscScanString()
383b3506946SBarry Smith .    Internal strings have arbitrary length and string copies are not checked that they fit into string space
384b3506946SBarry Smith -    Only works for PetscInt == int, PetscReal == double etc
385b3506946SBarry Smith 
386b3506946SBarry Smith 
387b3506946SBarry Smith */
388475446a1SBarry Smith PetscErrorCode PetscOptionsSAWsInput()
389b3506946SBarry Smith {
390b3506946SBarry Smith   PetscErrorCode ierr;
391b3506946SBarry Smith   PetscOptions   next     = PetscOptionsObject.next;
392d5649816SBarry Smith   static int     mancount = 0;
393b3506946SBarry Smith   char           options[16];
3947aab2a10SBarry Smith   PetscBool      changedmethod = PETSC_FALSE;
39588a9752cSBarry Smith   char           manname[16],textname[16];
3962657e9d9SBarry Smith   char           dir[1024];
397b3506946SBarry Smith 
398b3506946SBarry Smith   /* the next line is a bug, this will only work if all processors are here, the comm passed in is ignored!!! */
399b3506946SBarry Smith   sprintf(options,"Options_%d",count++);
400a297a907SKarl Rupp 
401e04113cfSBarry Smith   PetscOptionsObject.pprefix = PetscOptionsObject.prefix; /* SAWs will change this, so cannot pass prefix directly */
4021bc75a8dSBarry Smith 
403d50831c4SBarry Smith   ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s","_title");CHKERRQ(ierr);
4047781c08eSBarry Smith   PetscStackCallSAWs(SAWs_Register,(dir,&PetscOptionsObject.title,1,SAWs_READ,SAWs_STRING));
4057781c08eSBarry Smith   ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s","prefix");CHKERRQ(ierr);
4062657e9d9SBarry Smith   PetscStackCallSAWs(SAWs_Register,(dir,&PetscOptionsObject.pprefix,1,SAWs_READ,SAWs_STRING));
4072657e9d9SBarry Smith   PetscStackCallSAWs(SAWs_Register,("/PETSc/Options/ChangedMethod",&changedmethod,1,SAWs_WRITE,SAWs_BOOLEAN));
408b3506946SBarry Smith 
409b3506946SBarry Smith   while (next) {
410d50831c4SBarry Smith     sprintf(manname,"_man_%d",mancount);
4112657e9d9SBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",manname);CHKERRQ(ierr);
4127781c08eSBarry Smith     PetscStackCallSAWs(SAWs_Register,(dir,&next->man,1,SAWs_READ,SAWs_STRING));
413d50831c4SBarry Smith     sprintf(textname,"_text_%d",mancount++);
4147781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",textname);CHKERRQ(ierr);
4157781c08eSBarry Smith     PetscStackCallSAWs(SAWs_Register,(dir,&next->text,1,SAWs_READ,SAWs_STRING));
4169f32e415SBarry Smith 
417b3506946SBarry Smith     switch (next->type) {
418b3506946SBarry Smith     case OPTION_HEAD:
419b3506946SBarry Smith       break;
420b3506946SBarry Smith     case OPTION_INT_ARRAY:
4217781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4222657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_INT));
423b3506946SBarry Smith       break;
424b3506946SBarry Smith     case OPTION_REAL_ARRAY:
4257781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4262657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_DOUBLE));
427b3506946SBarry Smith       break;
428b3506946SBarry Smith     case OPTION_INT:
4297781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4302657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_INT));
431b3506946SBarry Smith       break;
432b3506946SBarry Smith     case OPTION_REAL:
4337781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4342657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_DOUBLE));
435b3506946SBarry Smith       break;
4367781c08eSBarry Smith     case OPTION_BOOL:
4377781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4382657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_BOOLEAN));
4391ae3d29cSBarry Smith       break;
4407781c08eSBarry Smith     case OPTION_BOOL_ARRAY:
4417781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4422657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_BOOLEAN));
44371f08665SBarry Smith       break;
444b3506946SBarry Smith     case OPTION_STRING:
4457781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4467781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
4471ae3d29cSBarry Smith       break;
4481ae3d29cSBarry Smith     case OPTION_STRING_ARRAY:
4497781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4502657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_STRING));
451b3506946SBarry Smith       break;
452a264d7a6SBarry Smith     case OPTION_FLIST:
453a264d7a6SBarry Smith       {
454a264d7a6SBarry Smith       PetscInt ntext;
4557781c08eSBarry Smith       ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4567781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
457a264d7a6SBarry Smith       ierr = PetscFunctionListGet(next->flist,(const char***)&next->edata,&ntext);CHKERRQ(ierr);
458a264d7a6SBarry Smith       PetscStackCallSAWs(SAWs_Set_Legal_Variable_Values,(dir,ntext,next->edata));
459a264d7a6SBarry Smith       }
460a264d7a6SBarry Smith       break;
4611ae3d29cSBarry Smith     case OPTION_ELIST:
462a264d7a6SBarry Smith       {
463a264d7a6SBarry Smith       PetscInt ntext = next->nlist;
4647781c08eSBarry Smith       ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4657781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
466785e854fSJed Brown       ierr = PetscMalloc1((ntext+1),&next->edata);CHKERRQ(ierr);
4671ae3d29cSBarry Smith       ierr = PetscMemcpy(next->edata,next->list,ntext*sizeof(char*));CHKERRQ(ierr);
468a264d7a6SBarry Smith       PetscStackCallSAWs(SAWs_Set_Legal_Variable_Values,(dir,ntext,next->edata));
469a264d7a6SBarry Smith       }
470a264d7a6SBarry Smith       break;
471b3506946SBarry Smith     default:
472b3506946SBarry Smith       break;
473b3506946SBarry Smith     }
474b3506946SBarry Smith     next = next->next;
475b3506946SBarry Smith   }
476b3506946SBarry Smith 
477b3506946SBarry Smith   /* wait until accessor has unlocked the memory */
4787aab2a10SBarry Smith   ierr = PetscSAWsBlock();CHKERRQ(ierr);
479b3506946SBarry Smith 
48088a9752cSBarry Smith   /* determine if any values have been set in GUI */
48188a9752cSBarry Smith   next = PetscOptionsObject.next;
48288a9752cSBarry Smith   while (next) {
48388a9752cSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
48488a9752cSBarry Smith     PetscStackCallSAWs(SAWs_Selected,(dir,&next->set));
48588a9752cSBarry Smith     next = next->next;
48688a9752cSBarry Smith   }
48788a9752cSBarry Smith 
488b3506946SBarry Smith   /* reset counter to -2; this updates the screen with the new options for the selected method */
489b3506946SBarry Smith   if (changedmethod) PetscOptionsPublishCount = -2;
490b3506946SBarry Smith 
4919a492a5cSBarry Smith   PetscStackCallSAWs(SAWs_Delete,("/PETSc/Options"));
492b3506946SBarry Smith   PetscFunctionReturn(0);
493b3506946SBarry Smith }
494b3506946SBarry Smith #endif
495b3506946SBarry Smith 
4966356e834SBarry Smith #undef __FUNCT__
49753acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEnd_Private"
49853acd3b1SBarry Smith PetscErrorCode PetscOptionsEnd_Private(void)
49953acd3b1SBarry Smith {
50053acd3b1SBarry Smith   PetscErrorCode ierr;
5016356e834SBarry Smith   PetscOptions   last;
5026356e834SBarry Smith   char           option[256],value[1024],tmp[32];
503330cf3c9SBarry Smith   size_t         j;
50453acd3b1SBarry Smith 
50553acd3b1SBarry Smith   PetscFunctionBegin;
506aee2cecaSBarry Smith   if (PetscOptionsObject.next) {
507b3506946SBarry Smith     if (!PetscOptionsPublishCount) {
508a264d7a6SBarry Smith #if defined(PETSC_HAVE_SAWS)
509475446a1SBarry Smith       ierr = PetscOptionsSAWsInput();CHKERRQ(ierr);
510b3506946SBarry Smith #else
51171f08665SBarry Smith       ierr = PetscOptionsGetFromTextInput();CHKERRQ(ierr);
512b3506946SBarry Smith #endif
513aee2cecaSBarry Smith     }
514aee2cecaSBarry Smith   }
5156356e834SBarry Smith 
516c31cb41cSBarry Smith   ierr = PetscFree(PetscOptionsObject.title);CHKERRQ(ierr);
517c31cb41cSBarry Smith   ierr = PetscFree(PetscOptionsObject.prefix);CHKERRQ(ierr);
5186356e834SBarry Smith 
5196356e834SBarry Smith   /* reset counter to -2; this updates the screen with the new options for the selected method */
5206356e834SBarry Smith   if (PetscOptionsObject.changedmethod) PetscOptionsPublishCount = -2;
52161b37b28SSatish Balay   /* reset alreadyprinted flag */
52261b37b28SSatish Balay   PetscOptionsObject.alreadyprinted = PETSC_FALSE;
5233194b578SJed Brown   if (PetscOptionsObject.object) PetscOptionsObject.object->optionsprinted = PETSC_TRUE;
5240298fd71SBarry Smith   PetscOptionsObject.object = NULL;
5256356e834SBarry Smith 
5266356e834SBarry Smith   while (PetscOptionsObject.next) {
5276356e834SBarry Smith     if (PetscOptionsObject.next->set) {
5286356e834SBarry Smith       if (PetscOptionsObject.prefix) {
5296356e834SBarry Smith         ierr = PetscStrcpy(option,"-");CHKERRQ(ierr);
5306356e834SBarry Smith         ierr = PetscStrcat(option,PetscOptionsObject.prefix);CHKERRQ(ierr);
5316356e834SBarry Smith         ierr = PetscStrcat(option,PetscOptionsObject.next->option+1);CHKERRQ(ierr);
5326356e834SBarry Smith       } else {
5336356e834SBarry Smith         ierr = PetscStrcpy(option,PetscOptionsObject.next->option);CHKERRQ(ierr);
5346356e834SBarry Smith       }
5356356e834SBarry Smith 
5366356e834SBarry Smith       switch (PetscOptionsObject.next->type) {
5376356e834SBarry Smith       case OPTION_HEAD:
5386356e834SBarry Smith         break;
539e26ddf31SBarry Smith       case OPTION_INT_ARRAY:
540e26ddf31SBarry Smith         sprintf(value,"%d",(int)((PetscInt*)PetscOptionsObject.next->data)[0]);
541e26ddf31SBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
542e26ddf31SBarry Smith           sprintf(tmp,"%d",(int)((PetscInt*)PetscOptionsObject.next->data)[j]);
543e26ddf31SBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
544e26ddf31SBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
545e26ddf31SBarry Smith         }
546e26ddf31SBarry Smith         break;
5476356e834SBarry Smith       case OPTION_INT:
5487a72a596SBarry Smith         sprintf(value,"%d",(int) *(PetscInt*)PetscOptionsObject.next->data);
5496356e834SBarry Smith         break;
5506356e834SBarry Smith       case OPTION_REAL:
5517a72a596SBarry Smith         sprintf(value,"%g",(double) *(PetscReal*)PetscOptionsObject.next->data);
5526356e834SBarry Smith         break;
5536356e834SBarry Smith       case OPTION_REAL_ARRAY:
5547a72a596SBarry Smith         sprintf(value,"%g",(double)((PetscReal*)PetscOptionsObject.next->data)[0]);
5556356e834SBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
5567a72a596SBarry Smith           sprintf(tmp,"%g",(double)((PetscReal*)PetscOptionsObject.next->data)[j]);
5576356e834SBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5586356e834SBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5596356e834SBarry Smith         }
5606356e834SBarry Smith         break;
5617781c08eSBarry Smith       case OPTION_BOOL:
56271f08665SBarry Smith         sprintf(value,"%d",*(int*)PetscOptionsObject.next->data);
5636356e834SBarry Smith         break;
5647781c08eSBarry Smith       case OPTION_BOOL_ARRAY:
565ace3abfcSBarry Smith         sprintf(value,"%d",(int)((PetscBool*)PetscOptionsObject.next->data)[0]);
5661ae3d29cSBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
567ace3abfcSBarry Smith           sprintf(tmp,"%d",(int)((PetscBool*)PetscOptionsObject.next->data)[j]);
5681ae3d29cSBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5691ae3d29cSBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5701ae3d29cSBarry Smith         }
5711ae3d29cSBarry Smith         break;
572a264d7a6SBarry Smith       case OPTION_FLIST:
5731ae3d29cSBarry Smith       case OPTION_ELIST:
5747781c08eSBarry Smith         ierr = PetscStrcpy(value,(char*)PetscOptionsObject.next->data);CHKERRQ(ierr);
5756356e834SBarry Smith         break;
5761ae3d29cSBarry Smith       case OPTION_STRING:
577475446a1SBarry Smith         ierr = PetscStrcpy(value,(char*)PetscOptionsObject.next->data);CHKERRQ(ierr);
5788da2146eSBarry Smith         break;
5791ae3d29cSBarry Smith       case OPTION_STRING_ARRAY:
5801ae3d29cSBarry Smith         sprintf(value,"%s",((char**)PetscOptionsObject.next->data)[0]);
5811ae3d29cSBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
5821ae3d29cSBarry Smith           sprintf(tmp,"%s",((char**)PetscOptionsObject.next->data)[j]);
5831ae3d29cSBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5841ae3d29cSBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5851ae3d29cSBarry Smith         }
5866356e834SBarry Smith         break;
5876356e834SBarry Smith       }
5886356e834SBarry Smith       ierr = PetscOptionsSetValue(option,value);CHKERRQ(ierr);
5896356e834SBarry Smith     }
590503cfb0cSBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->text);CHKERRQ(ierr);
591503cfb0cSBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->option);CHKERRQ(ierr);
5926356e834SBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->man);CHKERRQ(ierr);
59371f08665SBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->edata);CHKERRQ(ierr);
594c979a496SBarry Smith 
595c979a496SBarry Smith     if ((PetscOptionsObject.next->type == OPTION_STRING) || (PetscOptionsObject.next->type == OPTION_FLIST) || (PetscOptionsObject.next->type == OPTION_ELIST)){
59664facd6cSBarry Smith       /* must use system free since SAWs may have allocated it */
597c979a496SBarry Smith       free(PetscOptionsObject.next->data);
598c979a496SBarry Smith     } else {
5997781c08eSBarry Smith       ierr   = PetscFree(PetscOptionsObject.next->data);CHKERRQ(ierr);
600c979a496SBarry Smith     }
6017781c08eSBarry Smith 
6026356e834SBarry Smith     last                    = PetscOptionsObject.next;
6036356e834SBarry Smith     PetscOptionsObject.next = PetscOptionsObject.next->next;
6046356e834SBarry Smith     ierr                    = PetscFree(last);CHKERRQ(ierr);
6056356e834SBarry Smith   }
6066356e834SBarry Smith   PetscOptionsObject.next = 0;
60753acd3b1SBarry Smith   PetscFunctionReturn(0);
60853acd3b1SBarry Smith }
60953acd3b1SBarry Smith 
61053acd3b1SBarry Smith #undef __FUNCT__
61153acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEnum"
61253acd3b1SBarry Smith /*@C
61353acd3b1SBarry Smith    PetscOptionsEnum - Gets the enum value for a particular option in the database.
61453acd3b1SBarry Smith 
6153f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
61653acd3b1SBarry Smith 
61753acd3b1SBarry Smith    Input Parameters:
61853acd3b1SBarry Smith +  opt - option name
61953acd3b1SBarry Smith .  text - short string that describes the option
62053acd3b1SBarry Smith .  man - manual page with additional information on option
62153acd3b1SBarry Smith .  list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null
62253acd3b1SBarry Smith -  defaultv - the default (current) value
62353acd3b1SBarry Smith 
62453acd3b1SBarry Smith    Output Parameter:
62553acd3b1SBarry Smith +  value - the  value to return
626b32e0204SMatthew G Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
62753acd3b1SBarry Smith 
62853acd3b1SBarry Smith    Level: beginner
62953acd3b1SBarry Smith 
63053acd3b1SBarry Smith    Concepts: options database
63153acd3b1SBarry Smith 
63253acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
63353acd3b1SBarry Smith 
63453acd3b1SBarry Smith           list is usually something like PCASMTypes or some other predefined list of enum names
63553acd3b1SBarry Smith 
63653acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
637acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
638acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
63953acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
64053acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
641acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
642a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
64353acd3b1SBarry Smith @*/
6447087cfbeSBarry Smith PetscErrorCode  PetscOptionsEnum(const char opt[],const char text[],const char man[],const char *const *list,PetscEnum defaultv,PetscEnum *value,PetscBool  *set)
64553acd3b1SBarry Smith {
64653acd3b1SBarry Smith   PetscErrorCode ierr;
64753acd3b1SBarry Smith   PetscInt       ntext = 0;
648aa5bb8c0SSatish Balay   PetscInt       tval;
649ace3abfcSBarry Smith   PetscBool      tflg;
65053acd3b1SBarry Smith 
65153acd3b1SBarry Smith   PetscFunctionBegin;
65253acd3b1SBarry Smith   while (list[ntext++]) {
653e32f2f54SBarry Smith     if (ntext > 50) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument appears to be wrong or have more than 50 entries");
65453acd3b1SBarry Smith   }
655e32f2f54SBarry Smith   if (ntext < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument must have at least two entries: typename and type prefix");
65653acd3b1SBarry Smith   ntext -= 3;
657aa5bb8c0SSatish Balay   ierr   = PetscOptionsEList(opt,text,man,list,ntext,list[defaultv],&tval,&tflg);CHKERRQ(ierr);
658aa5bb8c0SSatish Balay   /* with PETSC_USE_64BIT_INDICES sizeof(PetscInt) != sizeof(PetscEnum) */
659aa5bb8c0SSatish Balay   if (tflg) *value = (PetscEnum)tval;
660aa5bb8c0SSatish Balay   if (set)  *set   = tflg;
66153acd3b1SBarry Smith   PetscFunctionReturn(0);
66253acd3b1SBarry Smith }
66353acd3b1SBarry Smith 
66453acd3b1SBarry Smith /* -------------------------------------------------------------------------------------------------------------*/
66553acd3b1SBarry Smith #undef __FUNCT__
66653acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsInt"
66753acd3b1SBarry Smith /*@C
66853acd3b1SBarry Smith    PetscOptionsInt - Gets the integer value for a particular option in the database.
66953acd3b1SBarry Smith 
6703f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
67153acd3b1SBarry Smith 
67253acd3b1SBarry Smith    Input Parameters:
67353acd3b1SBarry Smith +  opt - option name
67453acd3b1SBarry Smith .  text - short string that describes the option
67553acd3b1SBarry Smith .  man - manual page with additional information on option
676868c398cSBarry Smith -  defaultv - the default (current) value, if the user does not provide a value this is returned in value
67753acd3b1SBarry Smith 
67853acd3b1SBarry Smith    Output Parameter:
67953acd3b1SBarry Smith +  value - the integer value to return
68053acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
68153acd3b1SBarry Smith 
68253acd3b1SBarry Smith    Level: beginner
68353acd3b1SBarry Smith 
68453acd3b1SBarry Smith    Concepts: options database^has int
68553acd3b1SBarry Smith 
68653acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
68753acd3b1SBarry Smith 
68853acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
689acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
690acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
69153acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
69253acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
693acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
694a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
69553acd3b1SBarry Smith @*/
6967087cfbeSBarry Smith PetscErrorCode  PetscOptionsInt(const char opt[],const char text[],const char man[],PetscInt defaultv,PetscInt *value,PetscBool  *set)
69753acd3b1SBarry Smith {
69853acd3b1SBarry Smith   PetscErrorCode ierr;
6996356e834SBarry Smith   PetscOptions   amsopt;
70053acd3b1SBarry Smith 
70153acd3b1SBarry Smith   PetscFunctionBegin;
702b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
7036356e834SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_INT,&amsopt);CHKERRQ(ierr);
7046356e834SBarry Smith     ierr = PetscMalloc(sizeof(PetscInt),&amsopt->data);CHKERRQ(ierr);
705a297a907SKarl Rupp 
7066356e834SBarry Smith     *(PetscInt*)amsopt->data = defaultv;
707af6d86caSBarry Smith   }
70853acd3b1SBarry Smith   ierr = PetscOptionsGetInt(PetscOptionsObject.prefix,opt,value,set);CHKERRQ(ierr);
70961b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
7102aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%d>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,defaultv,text,ManSection(man));CHKERRQ(ierr);
71153acd3b1SBarry Smith   }
71253acd3b1SBarry Smith   PetscFunctionReturn(0);
71353acd3b1SBarry Smith }
71453acd3b1SBarry Smith 
71553acd3b1SBarry Smith #undef __FUNCT__
71653acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsString"
71753acd3b1SBarry Smith /*@C
71853acd3b1SBarry Smith    PetscOptionsString - Gets the string value for a particular option in the database.
71953acd3b1SBarry Smith 
7203f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
72153acd3b1SBarry Smith 
72253acd3b1SBarry Smith    Input Parameters:
72353acd3b1SBarry Smith +  opt - option name
72453acd3b1SBarry Smith .  text - short string that describes the option
72553acd3b1SBarry Smith .  man - manual page with additional information on option
726868c398cSBarry Smith .  defaultv - the default (current) value, if the user does not provide a value this is returned in value
727bcbf2dc5SJed Brown -  len - length of the result string including null terminator
72853acd3b1SBarry Smith 
72953acd3b1SBarry Smith    Output Parameter:
73053acd3b1SBarry Smith +  value - the value to return
73153acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
73253acd3b1SBarry Smith 
73353acd3b1SBarry Smith    Level: beginner
73453acd3b1SBarry Smith 
73553acd3b1SBarry Smith    Concepts: options database^has int
73653acd3b1SBarry Smith 
73753acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
73853acd3b1SBarry Smith 
7397fccdfe4SBarry 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).
7407fccdfe4SBarry Smith 
74153acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
742acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
743acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsReal(), PetscOptionsBool(),
74453acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
74553acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
746acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
747a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
74853acd3b1SBarry Smith @*/
7497087cfbeSBarry Smith PetscErrorCode  PetscOptionsString(const char opt[],const char text[],const char man[],const char defaultv[],char value[],size_t len,PetscBool  *set)
75053acd3b1SBarry Smith {
75153acd3b1SBarry Smith   PetscErrorCode ierr;
752aee2cecaSBarry Smith   PetscOptions   amsopt;
75353acd3b1SBarry Smith 
75453acd3b1SBarry Smith   PetscFunctionBegin;
755b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
756aee2cecaSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
75764facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
7585b02f95dSBarry Smith     ierr = PetscStrdup(defaultv ? defaultv : "",(char**)&amsopt->data);CHKERRQ(ierr);
759af6d86caSBarry Smith   }
76053acd3b1SBarry Smith   ierr = PetscOptionsGetString(PetscOptionsObject.prefix,opt,value,len,set);CHKERRQ(ierr);
76161b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
7622aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%s>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,defaultv,text,ManSection(man));CHKERRQ(ierr);
76353acd3b1SBarry Smith   }
76453acd3b1SBarry Smith   PetscFunctionReturn(0);
76553acd3b1SBarry Smith }
76653acd3b1SBarry Smith 
76753acd3b1SBarry Smith #undef __FUNCT__
76853acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsReal"
76953acd3b1SBarry Smith /*@C
77053acd3b1SBarry Smith    PetscOptionsReal - Gets the PetscReal value for a particular option in the database.
77153acd3b1SBarry Smith 
7723f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
77353acd3b1SBarry Smith 
77453acd3b1SBarry Smith    Input Parameters:
77553acd3b1SBarry Smith +  opt - option name
77653acd3b1SBarry Smith .  text - short string that describes the option
77753acd3b1SBarry Smith .  man - manual page with additional information on option
778868c398cSBarry Smith -  defaultv - the default (current) value, if the user does not provide a value this is returned in value
77953acd3b1SBarry Smith 
78053acd3b1SBarry Smith    Output Parameter:
78153acd3b1SBarry Smith +  value - the value to return
78253acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
78353acd3b1SBarry Smith 
78453acd3b1SBarry Smith    Level: beginner
78553acd3b1SBarry Smith 
78653acd3b1SBarry Smith    Concepts: options database^has int
78753acd3b1SBarry Smith 
78853acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
78953acd3b1SBarry Smith 
79053acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
791acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
792acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
79353acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
79453acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
795acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
796a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
79753acd3b1SBarry Smith @*/
7987087cfbeSBarry Smith PetscErrorCode  PetscOptionsReal(const char opt[],const char text[],const char man[],PetscReal defaultv,PetscReal *value,PetscBool  *set)
79953acd3b1SBarry Smith {
80053acd3b1SBarry Smith   PetscErrorCode ierr;
801538aa990SBarry Smith   PetscOptions   amsopt;
80253acd3b1SBarry Smith 
80353acd3b1SBarry Smith   PetscFunctionBegin;
804b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
805538aa990SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_REAL,&amsopt);CHKERRQ(ierr);
806538aa990SBarry Smith     ierr = PetscMalloc(sizeof(PetscReal),&amsopt->data);CHKERRQ(ierr);
807a297a907SKarl Rupp 
808538aa990SBarry Smith     *(PetscReal*)amsopt->data = defaultv;
809538aa990SBarry Smith   }
81053acd3b1SBarry Smith   ierr = PetscOptionsGetReal(PetscOptionsObject.prefix,opt,value,set);CHKERRQ(ierr);
81161b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
81257622a8eSBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%g>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,(double)defaultv,text,ManSection(man));CHKERRQ(ierr);
81353acd3b1SBarry Smith   }
81453acd3b1SBarry Smith   PetscFunctionReturn(0);
81553acd3b1SBarry Smith }
81653acd3b1SBarry Smith 
81753acd3b1SBarry Smith #undef __FUNCT__
81853acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsScalar"
81953acd3b1SBarry Smith /*@C
82053acd3b1SBarry Smith    PetscOptionsScalar - Gets the scalar value for a particular option in the database.
82153acd3b1SBarry Smith 
8223f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
82353acd3b1SBarry Smith 
82453acd3b1SBarry Smith    Input Parameters:
82553acd3b1SBarry Smith +  opt - option name
82653acd3b1SBarry Smith .  text - short string that describes the option
82753acd3b1SBarry Smith .  man - manual page with additional information on option
828868c398cSBarry Smith -  defaultv - the default (current) value, if the user does not provide a value this is returned in value
82953acd3b1SBarry Smith 
83053acd3b1SBarry Smith    Output Parameter:
83153acd3b1SBarry Smith +  value - the value to return
83253acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
83353acd3b1SBarry Smith 
83453acd3b1SBarry Smith    Level: beginner
83553acd3b1SBarry Smith 
83653acd3b1SBarry Smith    Concepts: options database^has int
83753acd3b1SBarry Smith 
83853acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
83953acd3b1SBarry Smith 
84053acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
841acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
842acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
84353acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
84453acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
845acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
846a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
84753acd3b1SBarry Smith @*/
8487087cfbeSBarry Smith PetscErrorCode  PetscOptionsScalar(const char opt[],const char text[],const char man[],PetscScalar defaultv,PetscScalar *value,PetscBool  *set)
84953acd3b1SBarry Smith {
85053acd3b1SBarry Smith   PetscErrorCode ierr;
85153acd3b1SBarry Smith 
85253acd3b1SBarry Smith   PetscFunctionBegin;
85353acd3b1SBarry Smith #if !defined(PETSC_USE_COMPLEX)
85453acd3b1SBarry Smith   ierr = PetscOptionsReal(opt,text,man,defaultv,value,set);CHKERRQ(ierr);
85553acd3b1SBarry Smith #else
85653acd3b1SBarry Smith   ierr = PetscOptionsGetScalar(PetscOptionsObject.prefix,opt,value,set);CHKERRQ(ierr);
85753acd3b1SBarry Smith #endif
85853acd3b1SBarry Smith   PetscFunctionReturn(0);
85953acd3b1SBarry Smith }
86053acd3b1SBarry Smith 
86153acd3b1SBarry Smith #undef __FUNCT__
86253acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsName"
86353acd3b1SBarry Smith /*@C
86490d69ab7SBarry 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
86590d69ab7SBarry Smith                       its value is set to false.
86653acd3b1SBarry Smith 
8673f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
86853acd3b1SBarry Smith 
86953acd3b1SBarry Smith    Input Parameters:
87053acd3b1SBarry Smith +  opt - option name
87153acd3b1SBarry Smith .  text - short string that describes the option
87253acd3b1SBarry Smith -  man - manual page with additional information on option
87353acd3b1SBarry Smith 
87453acd3b1SBarry Smith    Output Parameter:
87553acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
87653acd3b1SBarry Smith 
87753acd3b1SBarry Smith    Level: beginner
87853acd3b1SBarry Smith 
87953acd3b1SBarry Smith    Concepts: options database^has int
88053acd3b1SBarry Smith 
88153acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
88253acd3b1SBarry Smith 
88353acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
884acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
885acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
88653acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
88753acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
888acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
889a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
89053acd3b1SBarry Smith @*/
8917087cfbeSBarry Smith PetscErrorCode  PetscOptionsName(const char opt[],const char text[],const char man[],PetscBool  *flg)
89253acd3b1SBarry Smith {
89353acd3b1SBarry Smith   PetscErrorCode ierr;
8941ae3d29cSBarry Smith   PetscOptions   amsopt;
89553acd3b1SBarry Smith 
89653acd3b1SBarry Smith   PetscFunctionBegin;
8971ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
8987781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
899ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
900a297a907SKarl Rupp 
901ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
9021ae3d29cSBarry Smith   }
90353acd3b1SBarry Smith   ierr = PetscOptionsHasName(PetscOptionsObject.prefix,opt,flg);CHKERRQ(ierr);
90461b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
9052aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
90653acd3b1SBarry Smith   }
90753acd3b1SBarry Smith   PetscFunctionReturn(0);
90853acd3b1SBarry Smith }
90953acd3b1SBarry Smith 
91053acd3b1SBarry Smith #undef __FUNCT__
911a264d7a6SBarry Smith #define __FUNCT__ "PetscOptionsFList"
91253acd3b1SBarry Smith /*@C
913a264d7a6SBarry Smith      PetscOptionsFList - Puts a list of option values that a single one may be selected from
91453acd3b1SBarry Smith 
9153f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
91653acd3b1SBarry Smith 
91753acd3b1SBarry Smith    Input Parameters:
91853acd3b1SBarry Smith +  opt - option name
91953acd3b1SBarry Smith .  text - short string that describes the option
92053acd3b1SBarry Smith .  man - manual page with additional information on option
92153acd3b1SBarry Smith .  list - the possible choices
922868c398cSBarry Smith .  defaultv - the default (current) value, if the user does not provide a value this is returned in value
9233cc1e11dSBarry Smith -  len - the length of the character array value
92453acd3b1SBarry Smith 
92553acd3b1SBarry Smith    Output Parameter:
92653acd3b1SBarry Smith +  value - the value to return
92753acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
92853acd3b1SBarry Smith 
92953acd3b1SBarry Smith    Level: intermediate
93053acd3b1SBarry Smith 
93153acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
93253acd3b1SBarry Smith 
93353acd3b1SBarry Smith    See PetscOptionsEList() for when the choices are given in a string array
93453acd3b1SBarry Smith 
93553acd3b1SBarry Smith    To get a listing of all currently specified options,
93688c29154SBarry Smith     see PetscOptionsView() or PetscOptionsGetAll()
93753acd3b1SBarry Smith 
938eabe10d7SBarry Smith    Developer Note: This cannot check for invalid selection because of things like MATAIJ that are not included in the list
939eabe10d7SBarry Smith 
94053acd3b1SBarry Smith    Concepts: options database^list
94153acd3b1SBarry Smith 
94253acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
943acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
94453acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
94553acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
946acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
947a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsEnum()
94853acd3b1SBarry Smith @*/
949a264d7a6SBarry Smith PetscErrorCode  PetscOptionsFList(const char opt[],const char ltext[],const char man[],PetscFunctionList list,const char defaultv[],char value[],size_t len,PetscBool  *set)
95053acd3b1SBarry Smith {
95153acd3b1SBarry Smith   PetscErrorCode ierr;
9523cc1e11dSBarry Smith   PetscOptions   amsopt;
95353acd3b1SBarry Smith 
95453acd3b1SBarry Smith   PetscFunctionBegin;
955b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
956a264d7a6SBarry Smith     ierr = PetscOptionsCreate_Private(opt,ltext,man,OPTION_FLIST,&amsopt);CHKERRQ(ierr);
95764facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
9585b02f95dSBarry Smith     ierr = PetscStrdup(defaultv ? defaultv : "",(char**)&amsopt->data);CHKERRQ(ierr);
9593cc1e11dSBarry Smith     amsopt->flist = list;
9603cc1e11dSBarry Smith   }
96153acd3b1SBarry Smith   ierr = PetscOptionsGetString(PetscOptionsObject.prefix,opt,value,len,set);CHKERRQ(ierr);
96261b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
963140e18c1SBarry Smith     ierr = PetscFunctionListPrintTypes(PetscOptionsObject.comm,stdout,PetscOptionsObject.prefix,opt,ltext,man,list,defaultv);CHKERRQ(ierr);CHKERRQ(ierr);
96453acd3b1SBarry Smith   }
96553acd3b1SBarry Smith   PetscFunctionReturn(0);
96653acd3b1SBarry Smith }
96753acd3b1SBarry Smith 
96853acd3b1SBarry Smith #undef __FUNCT__
96953acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEList"
97053acd3b1SBarry Smith /*@C
97153acd3b1SBarry Smith      PetscOptionsEList - Puts a list of option values that a single one may be selected from
97253acd3b1SBarry Smith 
9733f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
97453acd3b1SBarry Smith 
97553acd3b1SBarry Smith    Input Parameters:
97653acd3b1SBarry Smith +  opt - option name
97753acd3b1SBarry Smith .  ltext - short string that describes the option
97853acd3b1SBarry Smith .  man - manual page with additional information on option
979a264d7a6SBarry Smith .  list - the possible choices (one of these must be selected, anything else is invalid)
98053acd3b1SBarry Smith .  ntext - number of choices
981868c398cSBarry Smith -  defaultv - the default (current) value, if the user does not provide a value the index of defaultv is returned
98253acd3b1SBarry Smith 
98353acd3b1SBarry Smith    Output Parameter:
98453acd3b1SBarry Smith +  value - the index of the value to return
98553acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
98653acd3b1SBarry Smith 
98753acd3b1SBarry Smith    Level: intermediate
98853acd3b1SBarry Smith 
98953acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
99053acd3b1SBarry Smith 
991a264d7a6SBarry Smith    See PetscOptionsFList() for when the choices are given in a PetscFunctionList()
99253acd3b1SBarry Smith 
99353acd3b1SBarry Smith    Concepts: options database^list
99453acd3b1SBarry Smith 
99553acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
996acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
99753acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
99853acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
999acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1000a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEnum()
100153acd3b1SBarry Smith @*/
10027087cfbeSBarry 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)
100353acd3b1SBarry Smith {
100453acd3b1SBarry Smith   PetscErrorCode ierr;
100553acd3b1SBarry Smith   PetscInt       i;
10061ae3d29cSBarry Smith   PetscOptions   amsopt;
100753acd3b1SBarry Smith 
100853acd3b1SBarry Smith   PetscFunctionBegin;
10091ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
1010d5649816SBarry Smith     ierr = PetscOptionsCreate_Private(opt,ltext,man,OPTION_ELIST,&amsopt);CHKERRQ(ierr);
101164facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
10125b02f95dSBarry Smith     ierr = PetscStrdup(defaultv ? defaultv : "",(char**)&amsopt->data);CHKERRQ(ierr);
10131ae3d29cSBarry Smith     amsopt->list  = list;
10141ae3d29cSBarry Smith     amsopt->nlist = ntext;
10151ae3d29cSBarry Smith   }
101653acd3b1SBarry Smith   ierr = PetscOptionsGetEList(PetscOptionsObject.prefix,opt,list,ntext,value,set);CHKERRQ(ierr);
101761b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
101853acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%s> (choose one of)",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,defaultv);CHKERRQ(ierr);
101953acd3b1SBarry Smith     for (i=0; i<ntext; i++) {
102053acd3b1SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm," %s",list[i]);CHKERRQ(ierr);
102153acd3b1SBarry Smith     }
10222aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm," (%s)\n",ManSection(man));CHKERRQ(ierr);
102353acd3b1SBarry Smith   }
102453acd3b1SBarry Smith   PetscFunctionReturn(0);
102553acd3b1SBarry Smith }
102653acd3b1SBarry Smith 
102753acd3b1SBarry Smith #undef __FUNCT__
1028acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolGroupBegin"
102953acd3b1SBarry Smith /*@C
1030acfcf0e5SJed Brown      PetscOptionsBoolGroupBegin - First in a series of logical queries on the options database for
1031d5649816SBarry Smith        which at most a single value can be true.
103253acd3b1SBarry Smith 
10333f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
103453acd3b1SBarry Smith 
103553acd3b1SBarry Smith    Input Parameters:
103653acd3b1SBarry Smith +  opt - option name
103753acd3b1SBarry Smith .  text - short string that describes the option
103853acd3b1SBarry Smith -  man - manual page with additional information on option
103953acd3b1SBarry Smith 
104053acd3b1SBarry Smith    Output Parameter:
104153acd3b1SBarry Smith .  flg - whether that option was set or not
104253acd3b1SBarry Smith 
104353acd3b1SBarry Smith    Level: intermediate
104453acd3b1SBarry Smith 
104553acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
104653acd3b1SBarry Smith 
1047acfcf0e5SJed Brown    Must be followed by 0 or more PetscOptionsBoolGroup()s and PetscOptionsBoolGroupEnd()
104853acd3b1SBarry Smith 
104953acd3b1SBarry Smith     Concepts: options database^logical group
105053acd3b1SBarry Smith 
105153acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1052acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
105353acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
105453acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1055acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1056a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
105753acd3b1SBarry Smith @*/
10587087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolGroupBegin(const char opt[],const char text[],const char man[],PetscBool  *flg)
105953acd3b1SBarry Smith {
106053acd3b1SBarry Smith   PetscErrorCode ierr;
10611ae3d29cSBarry Smith   PetscOptions   amsopt;
106253acd3b1SBarry Smith 
106353acd3b1SBarry Smith   PetscFunctionBegin;
10641ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
10657781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1066ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1067a297a907SKarl Rupp 
1068ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
10691ae3d29cSBarry Smith   }
107068b16fdaSBarry Smith   *flg = PETSC_FALSE;
10710298fd71SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,NULL);CHKERRQ(ierr);
107261b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
107353acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  Pick at most one of -------------\n");CHKERRQ(ierr);
10742aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"    -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
107553acd3b1SBarry Smith   }
107653acd3b1SBarry Smith   PetscFunctionReturn(0);
107753acd3b1SBarry Smith }
107853acd3b1SBarry Smith 
107953acd3b1SBarry Smith #undef __FUNCT__
1080acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolGroup"
108153acd3b1SBarry Smith /*@C
1082acfcf0e5SJed Brown      PetscOptionsBoolGroup - One in a series of logical queries on the options database for
1083d5649816SBarry Smith        which at most a single value can be true.
108453acd3b1SBarry Smith 
10853f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
108653acd3b1SBarry Smith 
108753acd3b1SBarry Smith    Input Parameters:
108853acd3b1SBarry Smith +  opt - option name
108953acd3b1SBarry Smith .  text - short string that describes the option
109053acd3b1SBarry Smith -  man - manual page with additional information on option
109153acd3b1SBarry Smith 
109253acd3b1SBarry Smith    Output Parameter:
109353acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
109453acd3b1SBarry Smith 
109553acd3b1SBarry Smith    Level: intermediate
109653acd3b1SBarry Smith 
109753acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
109853acd3b1SBarry Smith 
1099acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin() and preceded a PetscOptionsBoolGroupEnd()
110053acd3b1SBarry Smith 
110153acd3b1SBarry Smith     Concepts: options database^logical group
110253acd3b1SBarry Smith 
110353acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1104acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
110553acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
110653acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1107acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1108a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
110953acd3b1SBarry Smith @*/
11107087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolGroup(const char opt[],const char text[],const char man[],PetscBool  *flg)
111153acd3b1SBarry Smith {
111253acd3b1SBarry Smith   PetscErrorCode ierr;
11131ae3d29cSBarry Smith   PetscOptions   amsopt;
111453acd3b1SBarry Smith 
111553acd3b1SBarry Smith   PetscFunctionBegin;
11161ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
11177781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1118ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1119a297a907SKarl Rupp 
1120ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
11211ae3d29cSBarry Smith   }
112217326d04SJed Brown   *flg = PETSC_FALSE;
11230298fd71SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,NULL);CHKERRQ(ierr);
112461b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
11252aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"    -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
112653acd3b1SBarry Smith   }
112753acd3b1SBarry Smith   PetscFunctionReturn(0);
112853acd3b1SBarry Smith }
112953acd3b1SBarry Smith 
113053acd3b1SBarry Smith #undef __FUNCT__
1131acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolGroupEnd"
113253acd3b1SBarry Smith /*@C
1133acfcf0e5SJed Brown      PetscOptionsBoolGroupEnd - Last in a series of logical queries on the options database for
1134d5649816SBarry Smith        which at most a single value can be true.
113553acd3b1SBarry Smith 
11363f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
113753acd3b1SBarry Smith 
113853acd3b1SBarry Smith    Input Parameters:
113953acd3b1SBarry Smith +  opt - option name
114053acd3b1SBarry Smith .  text - short string that describes the option
114153acd3b1SBarry Smith -  man - manual page with additional information on option
114253acd3b1SBarry Smith 
114353acd3b1SBarry Smith    Output Parameter:
114453acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
114553acd3b1SBarry Smith 
114653acd3b1SBarry Smith    Level: intermediate
114753acd3b1SBarry Smith 
114853acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
114953acd3b1SBarry Smith 
1150acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin()
115153acd3b1SBarry Smith 
115253acd3b1SBarry Smith     Concepts: options database^logical group
115353acd3b1SBarry Smith 
115453acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1155acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
115653acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
115753acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1158acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1159a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
116053acd3b1SBarry Smith @*/
11617087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolGroupEnd(const char opt[],const char text[],const char man[],PetscBool  *flg)
116253acd3b1SBarry Smith {
116353acd3b1SBarry Smith   PetscErrorCode ierr;
11641ae3d29cSBarry Smith   PetscOptions   amsopt;
116553acd3b1SBarry Smith 
116653acd3b1SBarry Smith   PetscFunctionBegin;
11671ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
11687781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1169ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1170a297a907SKarl Rupp 
1171ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
11721ae3d29cSBarry Smith   }
117317326d04SJed Brown   *flg = PETSC_FALSE;
11740298fd71SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,NULL);CHKERRQ(ierr);
117561b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
11762aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"    -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
117753acd3b1SBarry Smith   }
117853acd3b1SBarry Smith   PetscFunctionReturn(0);
117953acd3b1SBarry Smith }
118053acd3b1SBarry Smith 
118153acd3b1SBarry Smith #undef __FUNCT__
1182acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBool"
118353acd3b1SBarry Smith /*@C
1184acfcf0e5SJed Brown    PetscOptionsBool - Determines if a particular option is in the database with a true or false
118553acd3b1SBarry Smith 
11863f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
118753acd3b1SBarry Smith 
118853acd3b1SBarry Smith    Input Parameters:
118953acd3b1SBarry Smith +  opt - option name
119053acd3b1SBarry Smith .  text - short string that describes the option
1191868c398cSBarry Smith .  man - manual page with additional information on option
1192868c398cSBarry Smith -  deflt - the default value, if the user does not set a value then this value is returned in flg
119353acd3b1SBarry Smith 
119453acd3b1SBarry Smith    Output Parameter:
119553acd3b1SBarry Smith .  flg - PETSC_TRUE or PETSC_FALSE
119653acd3b1SBarry Smith .  set - PETSC_TRUE if found, else PETSC_FALSE
119753acd3b1SBarry Smith 
119853acd3b1SBarry Smith    Level: beginner
119953acd3b1SBarry Smith 
120053acd3b1SBarry Smith    Concepts: options database^logical
120153acd3b1SBarry Smith 
120253acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
120353acd3b1SBarry Smith 
120453acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
1205acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
1206acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
120753acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
120853acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1209acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1210a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
121153acd3b1SBarry Smith @*/
12127087cfbeSBarry Smith PetscErrorCode  PetscOptionsBool(const char opt[],const char text[],const char man[],PetscBool deflt,PetscBool  *flg,PetscBool  *set)
121353acd3b1SBarry Smith {
121453acd3b1SBarry Smith   PetscErrorCode ierr;
1215ace3abfcSBarry Smith   PetscBool      iset;
1216aee2cecaSBarry Smith   PetscOptions   amsopt;
121753acd3b1SBarry Smith 
121853acd3b1SBarry Smith   PetscFunctionBegin;
1219b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
12207781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1221ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1222a297a907SKarl Rupp 
1223ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = deflt;
1224af6d86caSBarry Smith   }
1225acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,&iset);CHKERRQ(ierr);
122653acd3b1SBarry Smith   if (!iset) {
122753acd3b1SBarry Smith     if (flg) *flg = deflt;
122853acd3b1SBarry Smith   }
122953acd3b1SBarry Smith   if (set) *set = iset;
123061b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
1231ace3abfcSBarry Smith     const char *v = PetscBools[deflt];
12322aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s: <%s> %s (%s)\n",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,v,text,ManSection(man));CHKERRQ(ierr);
123353acd3b1SBarry Smith   }
123453acd3b1SBarry Smith   PetscFunctionReturn(0);
123553acd3b1SBarry Smith }
123653acd3b1SBarry Smith 
123753acd3b1SBarry Smith #undef __FUNCT__
123853acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsRealArray"
123953acd3b1SBarry Smith /*@C
124053acd3b1SBarry Smith    PetscOptionsRealArray - Gets an array of double values for a particular
124153acd3b1SBarry Smith    option in the database. The values must be separated with commas with
124253acd3b1SBarry Smith    no intervening spaces.
124353acd3b1SBarry Smith 
12443f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
124553acd3b1SBarry Smith 
124653acd3b1SBarry Smith    Input Parameters:
124753acd3b1SBarry Smith +  opt - the option one is seeking
124853acd3b1SBarry Smith .  text - short string describing option
124953acd3b1SBarry Smith .  man - manual page for option
125053acd3b1SBarry Smith -  nmax - maximum number of values
125153acd3b1SBarry Smith 
125253acd3b1SBarry Smith    Output Parameter:
125353acd3b1SBarry Smith +  value - location to copy values
125453acd3b1SBarry Smith .  nmax - actual number of values found
125553acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
125653acd3b1SBarry Smith 
125753acd3b1SBarry Smith    Level: beginner
125853acd3b1SBarry Smith 
125953acd3b1SBarry Smith    Notes:
126053acd3b1SBarry Smith    The user should pass in an array of doubles
126153acd3b1SBarry Smith 
126253acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
126353acd3b1SBarry Smith 
126453acd3b1SBarry Smith    Concepts: options database^array of strings
126553acd3b1SBarry Smith 
126653acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1267acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
126853acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
126953acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1270acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1271a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
127253acd3b1SBarry Smith @*/
12737087cfbeSBarry Smith PetscErrorCode  PetscOptionsRealArray(const char opt[],const char text[],const char man[],PetscReal value[],PetscInt *n,PetscBool  *set)
127453acd3b1SBarry Smith {
127553acd3b1SBarry Smith   PetscErrorCode ierr;
127653acd3b1SBarry Smith   PetscInt       i;
1277e26ddf31SBarry Smith   PetscOptions   amsopt;
127853acd3b1SBarry Smith 
127953acd3b1SBarry Smith   PetscFunctionBegin;
1280b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
1281e26ddf31SBarry Smith     PetscReal *vals;
1282e26ddf31SBarry Smith 
1283e26ddf31SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_REAL_ARRAY,&amsopt);CHKERRQ(ierr);
12846e02237eSPeter Brune     ierr = PetscMalloc1((*n),(PetscReal**)&amsopt->data);CHKERRQ(ierr);
1285e26ddf31SBarry Smith     vals = (PetscReal*)amsopt->data;
1286e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1287e26ddf31SBarry Smith     amsopt->arraylength = *n;
1288e26ddf31SBarry Smith   }
128953acd3b1SBarry Smith   ierr = PetscOptionsGetRealArray(PetscOptionsObject.prefix,opt,value,n,set);CHKERRQ(ierr);
129061b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
129157622a8eSBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%g",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,(double)value[0]);CHKERRQ(ierr);
129253acd3b1SBarry Smith     for (i=1; i<*n; i++) {
129357622a8eSBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,",%g",(double)value[i]);CHKERRQ(ierr);
129453acd3b1SBarry Smith     }
12952aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
129653acd3b1SBarry Smith   }
129753acd3b1SBarry Smith   PetscFunctionReturn(0);
129853acd3b1SBarry Smith }
129953acd3b1SBarry Smith 
130053acd3b1SBarry Smith 
130153acd3b1SBarry Smith #undef __FUNCT__
130253acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsIntArray"
130353acd3b1SBarry Smith /*@C
130453acd3b1SBarry Smith    PetscOptionsIntArray - Gets an array of integers for a particular
1305b32a342fSShri Abhyankar    option in the database.
130653acd3b1SBarry Smith 
13073f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
130853acd3b1SBarry Smith 
130953acd3b1SBarry Smith    Input Parameters:
131053acd3b1SBarry Smith +  opt - the option one is seeking
131153acd3b1SBarry Smith .  text - short string describing option
131253acd3b1SBarry Smith .  man - manual page for option
1313f8a50e2bSBarry Smith -  n - maximum number of values
131453acd3b1SBarry Smith 
131553acd3b1SBarry Smith    Output Parameter:
131653acd3b1SBarry Smith +  value - location to copy values
1317f8a50e2bSBarry Smith .  n - actual number of values found
131853acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
131953acd3b1SBarry Smith 
132053acd3b1SBarry Smith    Level: beginner
132153acd3b1SBarry Smith 
132253acd3b1SBarry Smith    Notes:
1323b32a342fSShri Abhyankar    The array can be passed as
1324b32a342fSShri Abhyankar    a comma seperated list:                                 0,1,2,3,4,5,6,7
13250fd488f5SShri Abhyankar    a range (start-end+1):                                  0-8
13260fd488f5SShri Abhyankar    a range with given increment (start-end+1:inc):         0-7:2
13270fd488f5SShri Abhyankar    a combination of values and ranges seperated by commas: 0,1-8,8-15:2
1328b32a342fSShri Abhyankar 
1329b32a342fSShri Abhyankar    There must be no intervening spaces between the values.
133053acd3b1SBarry Smith 
133153acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
133253acd3b1SBarry Smith 
1333b32a342fSShri Abhyankar    Concepts: options database^array of ints
133453acd3b1SBarry Smith 
133553acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1336acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
133753acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
133853acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1339acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1340a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsRealArray()
134153acd3b1SBarry Smith @*/
13427087cfbeSBarry Smith PetscErrorCode  PetscOptionsIntArray(const char opt[],const char text[],const char man[],PetscInt value[],PetscInt *n,PetscBool  *set)
134353acd3b1SBarry Smith {
134453acd3b1SBarry Smith   PetscErrorCode ierr;
134553acd3b1SBarry Smith   PetscInt       i;
1346e26ddf31SBarry Smith   PetscOptions   amsopt;
134753acd3b1SBarry Smith 
134853acd3b1SBarry Smith   PetscFunctionBegin;
1349b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
1350e26ddf31SBarry Smith     PetscInt *vals;
1351e26ddf31SBarry Smith 
1352e26ddf31SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_INT_ARRAY,&amsopt);CHKERRQ(ierr);
13536e02237eSPeter Brune     ierr = PetscMalloc1((*n),(PetscInt**)&amsopt->data);CHKERRQ(ierr);
1354e26ddf31SBarry Smith     vals = (PetscInt*)amsopt->data;
1355e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1356e26ddf31SBarry Smith     amsopt->arraylength = *n;
1357e26ddf31SBarry Smith   }
135853acd3b1SBarry Smith   ierr = PetscOptionsGetIntArray(PetscOptionsObject.prefix,opt,value,n,set);CHKERRQ(ierr);
135961b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
136053acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%d",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,value[0]);CHKERRQ(ierr);
136153acd3b1SBarry Smith     for (i=1; i<*n; i++) {
136253acd3b1SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,",%d",value[i]);CHKERRQ(ierr);
136353acd3b1SBarry Smith     }
13642aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
136553acd3b1SBarry Smith   }
136653acd3b1SBarry Smith   PetscFunctionReturn(0);
136753acd3b1SBarry Smith }
136853acd3b1SBarry Smith 
136953acd3b1SBarry Smith #undef __FUNCT__
137053acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsStringArray"
137153acd3b1SBarry Smith /*@C
137253acd3b1SBarry Smith    PetscOptionsStringArray - Gets an array of string values for a particular
137353acd3b1SBarry Smith    option in the database. The values must be separated with commas with
137453acd3b1SBarry Smith    no intervening spaces.
137553acd3b1SBarry Smith 
13763f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
137753acd3b1SBarry Smith 
137853acd3b1SBarry Smith    Input Parameters:
137953acd3b1SBarry Smith +  opt - the option one is seeking
138053acd3b1SBarry Smith .  text - short string describing option
138153acd3b1SBarry Smith .  man - manual page for option
138253acd3b1SBarry Smith -  nmax - maximum number of strings
138353acd3b1SBarry Smith 
138453acd3b1SBarry Smith    Output Parameter:
138553acd3b1SBarry Smith +  value - location to copy strings
138653acd3b1SBarry Smith .  nmax - actual number of strings found
138753acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
138853acd3b1SBarry Smith 
138953acd3b1SBarry Smith    Level: beginner
139053acd3b1SBarry Smith 
139153acd3b1SBarry Smith    Notes:
139253acd3b1SBarry Smith    The user should pass in an array of pointers to char, to hold all the
139353acd3b1SBarry Smith    strings returned by this function.
139453acd3b1SBarry Smith 
139553acd3b1SBarry Smith    The user is responsible for deallocating the strings that are
139653acd3b1SBarry Smith    returned. The Fortran interface for this routine is not supported.
139753acd3b1SBarry Smith 
139853acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
139953acd3b1SBarry Smith 
140053acd3b1SBarry Smith    Concepts: options database^array of strings
140153acd3b1SBarry Smith 
140253acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1403acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
140453acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
140553acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1406acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1407a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
140853acd3b1SBarry Smith @*/
14097087cfbeSBarry Smith PetscErrorCode  PetscOptionsStringArray(const char opt[],const char text[],const char man[],char *value[],PetscInt *nmax,PetscBool  *set)
141053acd3b1SBarry Smith {
141153acd3b1SBarry Smith   PetscErrorCode ierr;
14121ae3d29cSBarry Smith   PetscOptions   amsopt;
141353acd3b1SBarry Smith 
141453acd3b1SBarry Smith   PetscFunctionBegin;
14151ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
14161ae3d29cSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_STRING_ARRAY,&amsopt);CHKERRQ(ierr);
14176e02237eSPeter Brune     ierr = PetscMalloc1((*nmax),(char**)&amsopt->data);CHKERRQ(ierr);
1418a297a907SKarl Rupp 
14191ae3d29cSBarry Smith     amsopt->arraylength = *nmax;
14201ae3d29cSBarry Smith   }
142153acd3b1SBarry Smith   ierr = PetscOptionsGetStringArray(PetscOptionsObject.prefix,opt,value,nmax,set);CHKERRQ(ierr);
142261b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
14232aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <string1,string2,...>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
142453acd3b1SBarry Smith   }
142553acd3b1SBarry Smith   PetscFunctionReturn(0);
142653acd3b1SBarry Smith }
142753acd3b1SBarry Smith 
1428e2446a98SMatthew Knepley #undef __FUNCT__
1429acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolArray"
1430e2446a98SMatthew Knepley /*@C
1431acfcf0e5SJed Brown    PetscOptionsBoolArray - Gets an array of logical values (true or false) for a particular
1432e2446a98SMatthew Knepley    option in the database. The values must be separated with commas with
1433e2446a98SMatthew Knepley    no intervening spaces.
1434e2446a98SMatthew Knepley 
14353f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
1436e2446a98SMatthew Knepley 
1437e2446a98SMatthew Knepley    Input Parameters:
1438e2446a98SMatthew Knepley +  opt - the option one is seeking
1439e2446a98SMatthew Knepley .  text - short string describing option
1440e2446a98SMatthew Knepley .  man - manual page for option
1441e2446a98SMatthew Knepley -  nmax - maximum number of values
1442e2446a98SMatthew Knepley 
1443e2446a98SMatthew Knepley    Output Parameter:
1444e2446a98SMatthew Knepley +  value - location to copy values
1445e2446a98SMatthew Knepley .  nmax - actual number of values found
1446e2446a98SMatthew Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
1447e2446a98SMatthew Knepley 
1448e2446a98SMatthew Knepley    Level: beginner
1449e2446a98SMatthew Knepley 
1450e2446a98SMatthew Knepley    Notes:
1451e2446a98SMatthew Knepley    The user should pass in an array of doubles
1452e2446a98SMatthew Knepley 
1453e2446a98SMatthew Knepley    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
1454e2446a98SMatthew Knepley 
1455e2446a98SMatthew Knepley    Concepts: options database^array of strings
1456e2446a98SMatthew Knepley 
1457e2446a98SMatthew Knepley .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1458acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
1459e2446a98SMatthew Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1460e2446a98SMatthew Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1461acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1462a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
1463e2446a98SMatthew Knepley @*/
14647087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolArray(const char opt[],const char text[],const char man[],PetscBool value[],PetscInt *n,PetscBool *set)
1465e2446a98SMatthew Knepley {
1466e2446a98SMatthew Knepley   PetscErrorCode ierr;
1467e2446a98SMatthew Knepley   PetscInt       i;
14681ae3d29cSBarry Smith   PetscOptions   amsopt;
1469e2446a98SMatthew Knepley 
1470e2446a98SMatthew Knepley   PetscFunctionBegin;
14711ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
1472ace3abfcSBarry Smith     PetscBool *vals;
14731ae3d29cSBarry Smith 
14747781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL_ARRAY,&amsopt);CHKERRQ(ierr);
14756e02237eSPeter Brune     ierr = PetscMalloc1((*n),(PetscBool**)&amsopt->data);CHKERRQ(ierr);
1476ace3abfcSBarry Smith     vals = (PetscBool*)amsopt->data;
14771ae3d29cSBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
14781ae3d29cSBarry Smith     amsopt->arraylength = *n;
14791ae3d29cSBarry Smith   }
1480acfcf0e5SJed Brown   ierr = PetscOptionsGetBoolArray(PetscOptionsObject.prefix,opt,value,n,set);CHKERRQ(ierr);
1481e2446a98SMatthew Knepley   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
1482e2446a98SMatthew Knepley     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%d",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,value[0]);CHKERRQ(ierr);
1483e2446a98SMatthew Knepley     for (i=1; i<*n; i++) {
1484e2446a98SMatthew Knepley       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,",%d",value[i]);CHKERRQ(ierr);
1485e2446a98SMatthew Knepley     }
14862aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
1487e2446a98SMatthew Knepley   }
1488e2446a98SMatthew Knepley   PetscFunctionReturn(0);
1489e2446a98SMatthew Knepley }
1490e2446a98SMatthew Knepley 
14918cc676e6SMatthew G Knepley #undef __FUNCT__
14928cc676e6SMatthew G Knepley #define __FUNCT__ "PetscOptionsViewer"
14938cc676e6SMatthew G Knepley /*@C
1494d1da0b69SBarry Smith    PetscOptionsViewer - Gets a viewer appropriate for the type indicated by the user
14958cc676e6SMatthew G Knepley 
14968cc676e6SMatthew G Knepley    Logically Collective on the communicator passed in PetscOptionsBegin()
14978cc676e6SMatthew G Knepley 
14988cc676e6SMatthew G Knepley    Input Parameters:
14998cc676e6SMatthew G Knepley +  opt - option name
15008cc676e6SMatthew G Knepley .  text - short string that describes the option
15018cc676e6SMatthew G Knepley -  man - manual page with additional information on option
15028cc676e6SMatthew G Knepley 
15038cc676e6SMatthew G Knepley    Output Parameter:
15048cc676e6SMatthew G Knepley +  viewer - the viewer
15058cc676e6SMatthew G Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
15068cc676e6SMatthew G Knepley 
15078cc676e6SMatthew G Knepley    Level: beginner
15088cc676e6SMatthew G Knepley 
15098cc676e6SMatthew G Knepley    Concepts: options database^has int
15108cc676e6SMatthew G Knepley 
15118cc676e6SMatthew G Knepley    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
15128cc676e6SMatthew G Knepley 
1513d1da0b69SBarry Smith    See PetscOptionsGetVieweer() for the format of the supplied viewer and its options
15148cc676e6SMatthew G Knepley 
15158cc676e6SMatthew G Knepley .seealso: PetscOptionsGetViewer(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
15168cc676e6SMatthew G Knepley           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
15178cc676e6SMatthew G Knepley           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
15188cc676e6SMatthew G Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
15198cc676e6SMatthew G Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
15208cc676e6SMatthew G Knepley           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1521a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
15228cc676e6SMatthew G Knepley @*/
1523cffb1e40SBarry Smith PetscErrorCode  PetscOptionsViewer(const char opt[],const char text[],const char man[],PetscViewer *viewer,PetscViewerFormat *format,PetscBool  *set)
15248cc676e6SMatthew G Knepley {
15258cc676e6SMatthew G Knepley   PetscErrorCode ierr;
15268cc676e6SMatthew G Knepley   PetscOptions   amsopt;
15278cc676e6SMatthew G Knepley 
15288cc676e6SMatthew G Knepley   PetscFunctionBegin;
15298cc676e6SMatthew G Knepley   if (!PetscOptionsPublishCount) {
15308cc676e6SMatthew G Knepley     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
153164facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
15325b02f95dSBarry Smith     ierr = PetscStrdup("",(char**)&amsopt->data);CHKERRQ(ierr);
15338cc676e6SMatthew G Knepley   }
1534cffb1e40SBarry Smith   ierr = PetscOptionsGetViewer(PetscOptionsObject.comm,PetscOptionsObject.prefix,opt,viewer,format,set);CHKERRQ(ierr);
15358cc676e6SMatthew G Knepley   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
15368cc676e6SMatthew G Knepley     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%s>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,"",text,ManSection(man));CHKERRQ(ierr);
15378cc676e6SMatthew G Knepley   }
15388cc676e6SMatthew G Knepley   PetscFunctionReturn(0);
15398cc676e6SMatthew G Knepley }
15408cc676e6SMatthew G Knepley 
154153acd3b1SBarry Smith 
154253acd3b1SBarry Smith #undef __FUNCT__
154353acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsHead"
154453acd3b1SBarry Smith /*@C
1545b52f573bSBarry Smith      PetscOptionsHead - Puts a heading before listing any more published options. Used, for example,
154653acd3b1SBarry Smith             in KSPSetFromOptions_GMRES().
154753acd3b1SBarry Smith 
15483f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
154953acd3b1SBarry Smith 
155053acd3b1SBarry Smith    Input Parameter:
155153acd3b1SBarry Smith .   head - the heading text
155253acd3b1SBarry Smith 
155353acd3b1SBarry Smith 
155453acd3b1SBarry Smith    Level: intermediate
155553acd3b1SBarry Smith 
155653acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
155753acd3b1SBarry Smith 
1558b52f573bSBarry Smith           Can be followed by a call to PetscOptionsTail() in the same function.
155953acd3b1SBarry Smith 
156053acd3b1SBarry Smith    Concepts: options database^subheading
156153acd3b1SBarry Smith 
156253acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1563acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
156453acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
156553acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1566acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1567a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
156853acd3b1SBarry Smith @*/
15697087cfbeSBarry Smith PetscErrorCode  PetscOptionsHead(const char head[])
157053acd3b1SBarry Smith {
157153acd3b1SBarry Smith   PetscErrorCode ierr;
157253acd3b1SBarry Smith 
157353acd3b1SBarry Smith   PetscFunctionBegin;
157461b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
157553acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  %s\n",head);CHKERRQ(ierr);
157653acd3b1SBarry Smith   }
157753acd3b1SBarry Smith   PetscFunctionReturn(0);
157853acd3b1SBarry Smith }
157953acd3b1SBarry Smith 
158053acd3b1SBarry Smith 
158153acd3b1SBarry Smith 
158253acd3b1SBarry Smith 
158353acd3b1SBarry Smith 
158453acd3b1SBarry Smith 
1585