xref: /petsc/src/sys/objects/aoptions.c (revision d1da0b6977ee87f1a77143268e9e8db764d4e8b3) !
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 
1495b02f95dSBarry Smith /*
1505b02f95dSBarry Smith     This is needed because certain strings may be freed by SAWs, hence we cannot use PetscStrallocpy()
1515b02f95dSBarry Smith */
1525b02f95dSBarry Smith static PetscErrorCode  PetscStrdup(const char s[],char *t[])
1535b02f95dSBarry Smith {
1545b02f95dSBarry Smith   PetscErrorCode ierr;
1555b02f95dSBarry Smith   size_t         len;
1565b02f95dSBarry Smith   char           *tmp = 0;
1575b02f95dSBarry Smith 
1585b02f95dSBarry Smith   PetscFunctionBegin;
1595b02f95dSBarry Smith   if (s) {
1605b02f95dSBarry Smith     ierr = PetscStrlen(s,&len);CHKERRQ(ierr);
161ee48884fSBarry Smith     tmp = (char*) malloc((len+1)*sizeof(char*));
1625b02f95dSBarry Smith     if (!tmp) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"No memory to duplicate string");
1635b02f95dSBarry Smith     ierr = PetscStrcpy(tmp,s);CHKERRQ(ierr);
1645b02f95dSBarry Smith   }
1655b02f95dSBarry Smith   *t = tmp;
1665b02f95dSBarry Smith   PetscFunctionReturn(0);
1675b02f95dSBarry Smith }
1685b02f95dSBarry Smith 
1695b02f95dSBarry Smith 
170aee2cecaSBarry Smith #undef __FUNCT__
171aee2cecaSBarry Smith #define __FUNCT__ "PetscOptionsGetFromTextInput"
172aee2cecaSBarry Smith /*
1733cc1e11dSBarry Smith     PetscOptionsGetFromTextInput - Presents all the PETSc Options processed by the program so the user may change them at runtime
174aee2cecaSBarry Smith 
175aee2cecaSBarry Smith     Notes: this isn't really practical, it is just to demonstrate the principle
176aee2cecaSBarry Smith 
1777781c08eSBarry Smith     A carriage return indicates no change from the default; but this like -ksp_monitor <stdout>  the default is actually not stdout the default
1787781c08eSBarry Smith     is to do nothing so to get it to use stdout you need to type stdout. This is kind of bug?
1797781c08eSBarry Smith 
180aee2cecaSBarry Smith     Bugs:
1817781c08eSBarry Smith +    All processes must traverse through the exact same set of option queries due to the call to PetscScanString()
1823cc1e11dSBarry Smith .    Internal strings have arbitrary length and string copies are not checked that they fit into string space
183aee2cecaSBarry Smith -    Only works for PetscInt == int, PetscReal == double etc
184aee2cecaSBarry Smith 
1853cc1e11dSBarry Smith     Developer Notes: Normally the GUI that presents the options the user and retrieves the values would be running in a different
1863cc1e11dSBarry Smith      address space and communicating with the PETSc program
1873cc1e11dSBarry Smith 
188aee2cecaSBarry Smith */
189aee2cecaSBarry Smith PetscErrorCode PetscOptionsGetFromTextInput()
1906356e834SBarry Smith {
1916356e834SBarry Smith   PetscErrorCode ierr;
1926356e834SBarry Smith   PetscOptions   next = PetscOptionsObject.next;
1936356e834SBarry Smith   char           str[512];
1947781c08eSBarry Smith   PetscBool      bid;
195a4404d99SBarry Smith   PetscReal      ir,*valr;
196330cf3c9SBarry Smith   PetscInt       *vald;
197330cf3c9SBarry Smith   size_t         i;
1986356e834SBarry Smith 
199e26ddf31SBarry Smith   ierr = (*PetscPrintf)(PETSC_COMM_WORLD,"%s -------------------------------------------------\n",PetscOptionsObject.title);CHKERRQ(ierr);
2006356e834SBarry Smith   while (next) {
2016356e834SBarry Smith     switch (next->type) {
2026356e834SBarry Smith     case OPTION_HEAD:
2036356e834SBarry Smith       break;
204e26ddf31SBarry Smith     case OPTION_INT_ARRAY:
205e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",next->option+1);CHKERRQ(ierr);
206e26ddf31SBarry Smith       vald = (PetscInt*) next->data;
207e26ddf31SBarry Smith       for (i=0; i<next->arraylength; i++) {
208e26ddf31SBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"%d",vald[i]);CHKERRQ(ierr);
209e26ddf31SBarry Smith         if (i < next->arraylength-1) {
210e26ddf31SBarry Smith           ierr = PetscPrintf(PETSC_COMM_WORLD,",");CHKERRQ(ierr);
211e26ddf31SBarry Smith         }
212e26ddf31SBarry Smith       }
213e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,">: %s (%s) ",next->text,next->man);CHKERRQ(ierr);
214e26ddf31SBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
215e26ddf31SBarry Smith       if (str[0]) {
216e26ddf31SBarry Smith         PetscToken token;
217e26ddf31SBarry Smith         PetscInt   n=0,nmax = next->arraylength,*dvalue = (PetscInt*)next->data,start,end;
218e26ddf31SBarry Smith         size_t     len;
219e26ddf31SBarry Smith         char       *value;
220ace3abfcSBarry Smith         PetscBool  foundrange;
221e26ddf31SBarry Smith 
222e26ddf31SBarry Smith         next->set = PETSC_TRUE;
223e26ddf31SBarry Smith         value     = str;
224e26ddf31SBarry Smith         ierr      = PetscTokenCreate(value,',',&token);CHKERRQ(ierr);
225e26ddf31SBarry Smith         ierr      = PetscTokenFind(token,&value);CHKERRQ(ierr);
226e26ddf31SBarry Smith         while (n < nmax) {
227e26ddf31SBarry Smith           if (!value) break;
228e26ddf31SBarry Smith 
229e26ddf31SBarry Smith           /* look for form  d-D where d and D are integers */
230e26ddf31SBarry Smith           foundrange = PETSC_FALSE;
231e26ddf31SBarry Smith           ierr       = PetscStrlen(value,&len);CHKERRQ(ierr);
232e26ddf31SBarry Smith           if (value[0] == '-') i=2;
233e26ddf31SBarry Smith           else i=1;
234330cf3c9SBarry Smith           for (;i<len; i++) {
235e26ddf31SBarry Smith             if (value[i] == '-') {
236e32f2f54SBarry Smith               if (i == len-1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Error in %D-th array entry %s\n",n,value);
237e26ddf31SBarry Smith               value[i] = 0;
238cfbddea1SSatish Balay               ierr     = PetscOptionsStringToInt(value,&start);CHKERRQ(ierr);
239cfbddea1SSatish Balay               ierr     = PetscOptionsStringToInt(value+i+1,&end);CHKERRQ(ierr);
240e32f2f54SBarry 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);
241e32f2f54SBarry 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);
242e26ddf31SBarry Smith               for (; start<end; start++) {
243e26ddf31SBarry Smith                 *dvalue = start; dvalue++;n++;
244e26ddf31SBarry Smith               }
245e26ddf31SBarry Smith               foundrange = PETSC_TRUE;
246e26ddf31SBarry Smith               break;
247e26ddf31SBarry Smith             }
248e26ddf31SBarry Smith           }
249e26ddf31SBarry Smith           if (!foundrange) {
250cfbddea1SSatish Balay             ierr = PetscOptionsStringToInt(value,dvalue);CHKERRQ(ierr);
251e26ddf31SBarry Smith             dvalue++;
252e26ddf31SBarry Smith             n++;
253e26ddf31SBarry Smith           }
254e26ddf31SBarry Smith           ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
255e26ddf31SBarry Smith         }
2568c74ee41SBarry Smith         ierr = PetscTokenDestroy(&token);CHKERRQ(ierr);
257e26ddf31SBarry Smith       }
258e26ddf31SBarry Smith       break;
259e26ddf31SBarry Smith     case OPTION_REAL_ARRAY:
260e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",next->option+1);CHKERRQ(ierr);
261e26ddf31SBarry Smith       valr = (PetscReal*) next->data;
262e26ddf31SBarry Smith       for (i=0; i<next->arraylength; i++) {
263e26ddf31SBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"%g",valr[i]);CHKERRQ(ierr);
264e26ddf31SBarry Smith         if (i < next->arraylength-1) {
265e26ddf31SBarry Smith           ierr = PetscPrintf(PETSC_COMM_WORLD,",");CHKERRQ(ierr);
266e26ddf31SBarry Smith         }
267e26ddf31SBarry Smith       }
268e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,">: %s (%s) ",next->text,next->man);CHKERRQ(ierr);
269e26ddf31SBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
270e26ddf31SBarry Smith       if (str[0]) {
271e26ddf31SBarry Smith         PetscToken token;
272e26ddf31SBarry Smith         PetscInt   n = 0,nmax = next->arraylength;
273e26ddf31SBarry Smith         PetscReal  *dvalue = (PetscReal*)next->data;
274e26ddf31SBarry Smith         char       *value;
275e26ddf31SBarry Smith 
276e26ddf31SBarry Smith         next->set = PETSC_TRUE;
277e26ddf31SBarry Smith         value     = str;
278e26ddf31SBarry Smith         ierr      = PetscTokenCreate(value,',',&token);CHKERRQ(ierr);
279e26ddf31SBarry Smith         ierr      = PetscTokenFind(token,&value);CHKERRQ(ierr);
280e26ddf31SBarry Smith         while (n < nmax) {
281e26ddf31SBarry Smith           if (!value) break;
282cfbddea1SSatish Balay           ierr = PetscOptionsStringToReal(value,dvalue);CHKERRQ(ierr);
283e26ddf31SBarry Smith           dvalue++;
284e26ddf31SBarry Smith           n++;
285e26ddf31SBarry Smith           ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
286e26ddf31SBarry Smith         }
2878c74ee41SBarry Smith         ierr = PetscTokenDestroy(&token);CHKERRQ(ierr);
288e26ddf31SBarry Smith       }
289e26ddf31SBarry Smith       break;
2906356e834SBarry Smith     case OPTION_INT:
291e26ddf31SBarry 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);
2923fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
2933fc1eb6aSBarry Smith       if (str[0]) {
294d25d7f95SJed Brown #if defined(PETSC_SIZEOF_LONG_LONG)
295d25d7f95SJed Brown         long long lid;
296d25d7f95SJed Brown         sscanf(str,"%lld",&lid);
297d25d7f95SJed 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);
298c272547aSJed Brown #else
299d25d7f95SJed Brown         long  lid;
300d25d7f95SJed Brown         sscanf(str,"%ld",&lid);
301d25d7f95SJed 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);
302c272547aSJed Brown #endif
303a297a907SKarl Rupp 
304d25d7f95SJed Brown         next->set = PETSC_TRUE;
305d25d7f95SJed Brown         *((PetscInt*)next->data) = (PetscInt)lid;
306aee2cecaSBarry Smith       }
307aee2cecaSBarry Smith       break;
308aee2cecaSBarry Smith     case OPTION_REAL:
309e26ddf31SBarry 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);
3103fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3113fc1eb6aSBarry Smith       if (str[0]) {
312ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
313a4404d99SBarry Smith         sscanf(str,"%e",&ir);
314ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL_DOUBLE)
315aee2cecaSBarry Smith         sscanf(str,"%le",&ir);
316ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL___FLOAT128)
317d9822059SBarry Smith         ir = strtoflt128(str,0);
318d9822059SBarry Smith #else
319513dbe71SLisandro Dalcin         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Unknown scalar type");
320a4404d99SBarry Smith #endif
321aee2cecaSBarry Smith         next->set                 = PETSC_TRUE;
322aee2cecaSBarry Smith         *((PetscReal*)next->data) = ir;
323aee2cecaSBarry Smith       }
324aee2cecaSBarry Smith       break;
3257781c08eSBarry Smith     case OPTION_BOOL:
3267781c08eSBarry 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);
3277781c08eSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3287781c08eSBarry Smith       if (str[0]) {
3297781c08eSBarry Smith         ierr = PetscOptionsStringToBool(str,&bid);CHKERRQ(ierr);
3307781c08eSBarry Smith         next->set = PETSC_TRUE;
3317781c08eSBarry Smith         *((PetscBool*)next->data) = bid;
3327781c08eSBarry Smith       }
3337781c08eSBarry Smith       break;
334aee2cecaSBarry Smith     case OPTION_STRING:
335e26ddf31SBarry 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);
3363fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3373fc1eb6aSBarry Smith       if (str[0]) {
338aee2cecaSBarry Smith         next->set = PETSC_TRUE;
33964facd6cSBarry Smith         /* must use system malloc since SAWs may free this */
3405b02f95dSBarry Smith         ierr = PetscStrdup(str,(char**)&next->data);CHKERRQ(ierr);
3416356e834SBarry Smith       }
3426356e834SBarry Smith       break;
343a264d7a6SBarry Smith     case OPTION_FLIST:
344140e18c1SBarry Smith       ierr = PetscFunctionListPrintTypes(PETSC_COMM_WORLD,stdout,PetscOptionsObject.prefix,next->option,next->text,next->man,next->flist,(char*)next->data);CHKERRQ(ierr);
3453cc1e11dSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3463cc1e11dSBarry Smith       if (str[0]) {
3473cc1e11dSBarry Smith         PetscOptionsObject.changedmethod = PETSC_TRUE;
3483cc1e11dSBarry Smith         next->set = PETSC_TRUE;
34964facd6cSBarry Smith         /* must use system malloc since SAWs may free this */
3505b02f95dSBarry Smith         ierr = PetscStrdup(str,(char**)&next->data);CHKERRQ(ierr);
3513cc1e11dSBarry Smith       }
3523cc1e11dSBarry Smith       break;
353b432afdaSMatthew Knepley     default:
354b432afdaSMatthew Knepley       break;
3556356e834SBarry Smith     }
3566356e834SBarry Smith     next = next->next;
3576356e834SBarry Smith   }
3586356e834SBarry Smith   PetscFunctionReturn(0);
3596356e834SBarry Smith }
3606356e834SBarry Smith 
361e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
362e04113cfSBarry Smith #include <petscviewersaws.h>
363d5649816SBarry Smith 
364d5649816SBarry Smith static int count = 0;
365d5649816SBarry Smith 
366b3506946SBarry Smith #undef __FUNCT__
367e04113cfSBarry Smith #define __FUNCT__ "PetscOptionsSAWsDestroy"
368e04113cfSBarry Smith PetscErrorCode PetscOptionsSAWsDestroy(void)
369d5649816SBarry Smith {
3702657e9d9SBarry Smith   PetscFunctionBegin;
371d5649816SBarry Smith   PetscFunctionReturn(0);
372d5649816SBarry Smith }
373d5649816SBarry Smith 
374d5649816SBarry Smith #undef __FUNCT__
3757781c08eSBarry Smith #define __FUNCT__ "PetscOptionsSAWsInput"
376b3506946SBarry Smith /*
3777781c08eSBarry Smith     PetscOptionsSAWsInput - Presents all the PETSc Options processed by the program so the user may change them at runtime using the SAWs
378b3506946SBarry Smith 
379b3506946SBarry Smith     Bugs:
380b3506946SBarry Smith +    All processes must traverse through the exact same set of option queries do to the call to PetscScanString()
381b3506946SBarry Smith .    Internal strings have arbitrary length and string copies are not checked that they fit into string space
382b3506946SBarry Smith -    Only works for PetscInt == int, PetscReal == double etc
383b3506946SBarry Smith 
384b3506946SBarry Smith 
385b3506946SBarry Smith */
386475446a1SBarry Smith PetscErrorCode PetscOptionsSAWsInput()
387b3506946SBarry Smith {
388b3506946SBarry Smith   PetscErrorCode ierr;
389b3506946SBarry Smith   PetscOptions   next     = PetscOptionsObject.next;
390d5649816SBarry Smith   static int     mancount = 0;
391b3506946SBarry Smith   char           options[16];
3927aab2a10SBarry Smith   PetscBool      changedmethod = PETSC_FALSE;
39388a9752cSBarry Smith   char           manname[16],textname[16];
3942657e9d9SBarry Smith   char           dir[1024];
395b3506946SBarry Smith 
396b3506946SBarry Smith   /* the next line is a bug, this will only work if all processors are here, the comm passed in is ignored!!! */
397b3506946SBarry Smith   sprintf(options,"Options_%d",count++);
398a297a907SKarl Rupp 
399e04113cfSBarry Smith   PetscOptionsObject.pprefix = PetscOptionsObject.prefix; /* SAWs will change this, so cannot pass prefix directly */
4001bc75a8dSBarry Smith 
401d50831c4SBarry Smith   ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s","_title");CHKERRQ(ierr);
4027781c08eSBarry Smith   PetscStackCallSAWs(SAWs_Register,(dir,&PetscOptionsObject.title,1,SAWs_READ,SAWs_STRING));
4037781c08eSBarry Smith   ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s","prefix");CHKERRQ(ierr);
4042657e9d9SBarry Smith   PetscStackCallSAWs(SAWs_Register,(dir,&PetscOptionsObject.pprefix,1,SAWs_READ,SAWs_STRING));
4052657e9d9SBarry Smith   PetscStackCallSAWs(SAWs_Register,("/PETSc/Options/ChangedMethod",&changedmethod,1,SAWs_WRITE,SAWs_BOOLEAN));
406b3506946SBarry Smith 
407b3506946SBarry Smith   while (next) {
408d50831c4SBarry Smith     sprintf(manname,"_man_%d",mancount);
4092657e9d9SBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",manname);CHKERRQ(ierr);
4107781c08eSBarry Smith     PetscStackCallSAWs(SAWs_Register,(dir,&next->man,1,SAWs_READ,SAWs_STRING));
411d50831c4SBarry Smith     sprintf(textname,"_text_%d",mancount++);
4127781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",textname);CHKERRQ(ierr);
4137781c08eSBarry Smith     PetscStackCallSAWs(SAWs_Register,(dir,&next->text,1,SAWs_READ,SAWs_STRING));
4149f32e415SBarry Smith 
415b3506946SBarry Smith     switch (next->type) {
416b3506946SBarry Smith     case OPTION_HEAD:
417b3506946SBarry Smith       break;
418b3506946SBarry Smith     case OPTION_INT_ARRAY:
4197781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4202657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_INT));
421b3506946SBarry Smith       break;
422b3506946SBarry Smith     case OPTION_REAL_ARRAY:
4237781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4242657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_DOUBLE));
425b3506946SBarry Smith       break;
426b3506946SBarry Smith     case OPTION_INT:
4277781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4282657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_INT));
429b3506946SBarry Smith       break;
430b3506946SBarry Smith     case OPTION_REAL:
4317781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4322657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_DOUBLE));
433b3506946SBarry Smith       break;
4347781c08eSBarry Smith     case OPTION_BOOL:
4357781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4362657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_BOOLEAN));
4371ae3d29cSBarry Smith       break;
4387781c08eSBarry Smith     case OPTION_BOOL_ARRAY:
4397781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4402657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_BOOLEAN));
44171f08665SBarry Smith       break;
442b3506946SBarry Smith     case OPTION_STRING:
4437781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4447781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
4451ae3d29cSBarry Smith       break;
4461ae3d29cSBarry Smith     case OPTION_STRING_ARRAY:
4477781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4482657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_STRING));
449b3506946SBarry Smith       break;
450a264d7a6SBarry Smith     case OPTION_FLIST:
451a264d7a6SBarry Smith       {
452a264d7a6SBarry Smith       PetscInt ntext;
4537781c08eSBarry Smith       ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4547781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
455a264d7a6SBarry Smith       ierr = PetscFunctionListGet(next->flist,(const char***)&next->edata,&ntext);CHKERRQ(ierr);
456a264d7a6SBarry Smith       PetscStackCallSAWs(SAWs_Set_Legal_Variable_Values,(dir,ntext,next->edata));
457a264d7a6SBarry Smith       }
458a264d7a6SBarry Smith       break;
4591ae3d29cSBarry Smith     case OPTION_ELIST:
460a264d7a6SBarry Smith       {
461a264d7a6SBarry Smith       PetscInt ntext = next->nlist;
4627781c08eSBarry Smith       ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4637781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
464785e854fSJed Brown       ierr = PetscMalloc1((ntext+1),&next->edata);CHKERRQ(ierr);
4651ae3d29cSBarry Smith       ierr = PetscMemcpy(next->edata,next->list,ntext*sizeof(char*));CHKERRQ(ierr);
466a264d7a6SBarry Smith       PetscStackCallSAWs(SAWs_Set_Legal_Variable_Values,(dir,ntext,next->edata));
467a264d7a6SBarry Smith       }
468a264d7a6SBarry Smith       break;
469b3506946SBarry Smith     default:
470b3506946SBarry Smith       break;
471b3506946SBarry Smith     }
472b3506946SBarry Smith     next = next->next;
473b3506946SBarry Smith   }
474b3506946SBarry Smith 
475b3506946SBarry Smith   /* wait until accessor has unlocked the memory */
4767aab2a10SBarry Smith   ierr = PetscSAWsBlock();CHKERRQ(ierr);
477b3506946SBarry Smith 
47888a9752cSBarry Smith   /* determine if any values have been set in GUI */
47988a9752cSBarry Smith   next = PetscOptionsObject.next;
48088a9752cSBarry Smith   while (next) {
48188a9752cSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
48288a9752cSBarry Smith     PetscStackCallSAWs(SAWs_Selected,(dir,&next->set));
48388a9752cSBarry Smith     next = next->next;
48488a9752cSBarry Smith   }
48588a9752cSBarry Smith 
486b3506946SBarry Smith   /* reset counter to -2; this updates the screen with the new options for the selected method */
487b3506946SBarry Smith   if (changedmethod) PetscOptionsPublishCount = -2;
488b3506946SBarry Smith 
4899a492a5cSBarry Smith   PetscStackCallSAWs(SAWs_Delete,("/PETSc/Options"));
490b3506946SBarry Smith   PetscFunctionReturn(0);
491b3506946SBarry Smith }
492b3506946SBarry Smith #endif
493b3506946SBarry Smith 
4946356e834SBarry Smith #undef __FUNCT__
49553acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEnd_Private"
49653acd3b1SBarry Smith PetscErrorCode PetscOptionsEnd_Private(void)
49753acd3b1SBarry Smith {
49853acd3b1SBarry Smith   PetscErrorCode ierr;
4996356e834SBarry Smith   PetscOptions   last;
5006356e834SBarry Smith   char           option[256],value[1024],tmp[32];
501330cf3c9SBarry Smith   size_t         j;
50253acd3b1SBarry Smith 
50353acd3b1SBarry Smith   PetscFunctionBegin;
504aee2cecaSBarry Smith   if (PetscOptionsObject.next) {
505b3506946SBarry Smith     if (!PetscOptionsPublishCount) {
506a264d7a6SBarry Smith #if defined(PETSC_HAVE_SAWS)
507475446a1SBarry Smith       ierr = PetscOptionsSAWsInput();CHKERRQ(ierr);
508b3506946SBarry Smith #else
50971f08665SBarry Smith       ierr = PetscOptionsGetFromTextInput();CHKERRQ(ierr);
510b3506946SBarry Smith #endif
511aee2cecaSBarry Smith     }
512aee2cecaSBarry Smith   }
5136356e834SBarry Smith 
514c31cb41cSBarry Smith   ierr = PetscFree(PetscOptionsObject.title);CHKERRQ(ierr);
515c31cb41cSBarry Smith   ierr = PetscFree(PetscOptionsObject.prefix);CHKERRQ(ierr);
5166356e834SBarry Smith 
5176356e834SBarry Smith   /* reset counter to -2; this updates the screen with the new options for the selected method */
5186356e834SBarry Smith   if (PetscOptionsObject.changedmethod) PetscOptionsPublishCount = -2;
51961b37b28SSatish Balay   /* reset alreadyprinted flag */
52061b37b28SSatish Balay   PetscOptionsObject.alreadyprinted = PETSC_FALSE;
5213194b578SJed Brown   if (PetscOptionsObject.object) PetscOptionsObject.object->optionsprinted = PETSC_TRUE;
5220298fd71SBarry Smith   PetscOptionsObject.object = NULL;
5236356e834SBarry Smith 
5246356e834SBarry Smith   while (PetscOptionsObject.next) {
5256356e834SBarry Smith     if (PetscOptionsObject.next->set) {
5266356e834SBarry Smith       if (PetscOptionsObject.prefix) {
5276356e834SBarry Smith         ierr = PetscStrcpy(option,"-");CHKERRQ(ierr);
5286356e834SBarry Smith         ierr = PetscStrcat(option,PetscOptionsObject.prefix);CHKERRQ(ierr);
5296356e834SBarry Smith         ierr = PetscStrcat(option,PetscOptionsObject.next->option+1);CHKERRQ(ierr);
5306356e834SBarry Smith       } else {
5316356e834SBarry Smith         ierr = PetscStrcpy(option,PetscOptionsObject.next->option);CHKERRQ(ierr);
5326356e834SBarry Smith       }
5336356e834SBarry Smith 
5346356e834SBarry Smith       switch (PetscOptionsObject.next->type) {
5356356e834SBarry Smith       case OPTION_HEAD:
5366356e834SBarry Smith         break;
537e26ddf31SBarry Smith       case OPTION_INT_ARRAY:
538e26ddf31SBarry Smith         sprintf(value,"%d",(int)((PetscInt*)PetscOptionsObject.next->data)[0]);
539e26ddf31SBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
540e26ddf31SBarry Smith           sprintf(tmp,"%d",(int)((PetscInt*)PetscOptionsObject.next->data)[j]);
541e26ddf31SBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
542e26ddf31SBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
543e26ddf31SBarry Smith         }
544e26ddf31SBarry Smith         break;
5456356e834SBarry Smith       case OPTION_INT:
5467a72a596SBarry Smith         sprintf(value,"%d",(int) *(PetscInt*)PetscOptionsObject.next->data);
5476356e834SBarry Smith         break;
5486356e834SBarry Smith       case OPTION_REAL:
5497a72a596SBarry Smith         sprintf(value,"%g",(double) *(PetscReal*)PetscOptionsObject.next->data);
5506356e834SBarry Smith         break;
5516356e834SBarry Smith       case OPTION_REAL_ARRAY:
5527a72a596SBarry Smith         sprintf(value,"%g",(double)((PetscReal*)PetscOptionsObject.next->data)[0]);
5536356e834SBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
5547a72a596SBarry Smith           sprintf(tmp,"%g",(double)((PetscReal*)PetscOptionsObject.next->data)[j]);
5556356e834SBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5566356e834SBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5576356e834SBarry Smith         }
5586356e834SBarry Smith         break;
5597781c08eSBarry Smith       case OPTION_BOOL:
56071f08665SBarry Smith         sprintf(value,"%d",*(int*)PetscOptionsObject.next->data);
5616356e834SBarry Smith         break;
5627781c08eSBarry Smith       case OPTION_BOOL_ARRAY:
563ace3abfcSBarry Smith         sprintf(value,"%d",(int)((PetscBool*)PetscOptionsObject.next->data)[0]);
5641ae3d29cSBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
565ace3abfcSBarry Smith           sprintf(tmp,"%d",(int)((PetscBool*)PetscOptionsObject.next->data)[j]);
5661ae3d29cSBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5671ae3d29cSBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5681ae3d29cSBarry Smith         }
5691ae3d29cSBarry Smith         break;
570a264d7a6SBarry Smith       case OPTION_FLIST:
5711ae3d29cSBarry Smith       case OPTION_ELIST:
5727781c08eSBarry Smith         ierr = PetscStrcpy(value,(char*)PetscOptionsObject.next->data);CHKERRQ(ierr);
5736356e834SBarry Smith         break;
5741ae3d29cSBarry Smith       case OPTION_STRING:
575475446a1SBarry Smith         ierr = PetscStrcpy(value,(char*)PetscOptionsObject.next->data);CHKERRQ(ierr);
5768da2146eSBarry Smith         break;
5771ae3d29cSBarry Smith       case OPTION_STRING_ARRAY:
5781ae3d29cSBarry Smith         sprintf(value,"%s",((char**)PetscOptionsObject.next->data)[0]);
5791ae3d29cSBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
5801ae3d29cSBarry Smith           sprintf(tmp,"%s",((char**)PetscOptionsObject.next->data)[j]);
5811ae3d29cSBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5821ae3d29cSBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5831ae3d29cSBarry Smith         }
5846356e834SBarry Smith         break;
5856356e834SBarry Smith       }
5866356e834SBarry Smith       ierr = PetscOptionsSetValue(option,value);CHKERRQ(ierr);
5876356e834SBarry Smith     }
588503cfb0cSBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->text);CHKERRQ(ierr);
589503cfb0cSBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->option);CHKERRQ(ierr);
5906356e834SBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->man);CHKERRQ(ierr);
59171f08665SBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->edata);CHKERRQ(ierr);
592c979a496SBarry Smith 
593c979a496SBarry Smith     if ((PetscOptionsObject.next->type == OPTION_STRING) || (PetscOptionsObject.next->type == OPTION_FLIST) || (PetscOptionsObject.next->type == OPTION_ELIST)){
59464facd6cSBarry Smith       /* must use system free since SAWs may have allocated it */
595c979a496SBarry Smith       free(PetscOptionsObject.next->data);
596c979a496SBarry Smith     } else {
5977781c08eSBarry Smith       ierr   = PetscFree(PetscOptionsObject.next->data);CHKERRQ(ierr);
598c979a496SBarry Smith     }
5997781c08eSBarry Smith 
6006356e834SBarry Smith     last                    = PetscOptionsObject.next;
6016356e834SBarry Smith     PetscOptionsObject.next = PetscOptionsObject.next->next;
6026356e834SBarry Smith     ierr                    = PetscFree(last);CHKERRQ(ierr);
6036356e834SBarry Smith   }
6046356e834SBarry Smith   PetscOptionsObject.next = 0;
60553acd3b1SBarry Smith   PetscFunctionReturn(0);
60653acd3b1SBarry Smith }
60753acd3b1SBarry Smith 
60853acd3b1SBarry Smith #undef __FUNCT__
60953acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEnum"
61053acd3b1SBarry Smith /*@C
61153acd3b1SBarry Smith    PetscOptionsEnum - Gets the enum value for a particular option in the database.
61253acd3b1SBarry Smith 
6133f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
61453acd3b1SBarry Smith 
61553acd3b1SBarry Smith    Input Parameters:
61653acd3b1SBarry Smith +  opt - option name
61753acd3b1SBarry Smith .  text - short string that describes the option
61853acd3b1SBarry Smith .  man - manual page with additional information on option
61953acd3b1SBarry Smith .  list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null
62053acd3b1SBarry Smith -  defaultv - the default (current) value
62153acd3b1SBarry Smith 
62253acd3b1SBarry Smith    Output Parameter:
62353acd3b1SBarry Smith +  value - the  value to return
624b32e0204SMatthew G Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
62553acd3b1SBarry Smith 
62653acd3b1SBarry Smith    Level: beginner
62753acd3b1SBarry Smith 
62853acd3b1SBarry Smith    Concepts: options database
62953acd3b1SBarry Smith 
63053acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
63153acd3b1SBarry Smith 
63253acd3b1SBarry Smith           list is usually something like PCASMTypes or some other predefined list of enum names
63353acd3b1SBarry Smith 
63453acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
635acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
636acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
63753acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
63853acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
639acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
640a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
64153acd3b1SBarry Smith @*/
6427087cfbeSBarry Smith PetscErrorCode  PetscOptionsEnum(const char opt[],const char text[],const char man[],const char *const *list,PetscEnum defaultv,PetscEnum *value,PetscBool  *set)
64353acd3b1SBarry Smith {
64453acd3b1SBarry Smith   PetscErrorCode ierr;
64553acd3b1SBarry Smith   PetscInt       ntext = 0;
646aa5bb8c0SSatish Balay   PetscInt       tval;
647ace3abfcSBarry Smith   PetscBool      tflg;
64853acd3b1SBarry Smith 
64953acd3b1SBarry Smith   PetscFunctionBegin;
65053acd3b1SBarry Smith   while (list[ntext++]) {
651e32f2f54SBarry Smith     if (ntext > 50) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument appears to be wrong or have more than 50 entries");
65253acd3b1SBarry Smith   }
653e32f2f54SBarry Smith   if (ntext < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument must have at least two entries: typename and type prefix");
65453acd3b1SBarry Smith   ntext -= 3;
655aa5bb8c0SSatish Balay   ierr   = PetscOptionsEList(opt,text,man,list,ntext,list[defaultv],&tval,&tflg);CHKERRQ(ierr);
656aa5bb8c0SSatish Balay   /* with PETSC_USE_64BIT_INDICES sizeof(PetscInt) != sizeof(PetscEnum) */
657aa5bb8c0SSatish Balay   if (tflg) *value = (PetscEnum)tval;
658aa5bb8c0SSatish Balay   if (set)  *set   = tflg;
65953acd3b1SBarry Smith   PetscFunctionReturn(0);
66053acd3b1SBarry Smith }
66153acd3b1SBarry Smith 
66253acd3b1SBarry Smith /* -------------------------------------------------------------------------------------------------------------*/
66353acd3b1SBarry Smith #undef __FUNCT__
66453acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsInt"
66553acd3b1SBarry Smith /*@C
66653acd3b1SBarry Smith    PetscOptionsInt - Gets the integer value for a particular option in the database.
66753acd3b1SBarry Smith 
6683f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
66953acd3b1SBarry Smith 
67053acd3b1SBarry Smith    Input Parameters:
67153acd3b1SBarry Smith +  opt - option name
67253acd3b1SBarry Smith .  text - short string that describes the option
67353acd3b1SBarry Smith .  man - manual page with additional information on option
67453acd3b1SBarry Smith -  defaultv - the default (current) value
67553acd3b1SBarry Smith 
67653acd3b1SBarry Smith    Output Parameter:
67753acd3b1SBarry Smith +  value - the integer value to return
67853acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
67953acd3b1SBarry Smith 
68053acd3b1SBarry Smith    Level: beginner
68153acd3b1SBarry Smith 
68253acd3b1SBarry Smith    Concepts: options database^has int
68353acd3b1SBarry Smith 
68453acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
68553acd3b1SBarry Smith 
68653acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
687acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
688acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
68953acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
69053acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
691acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
692a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
69353acd3b1SBarry Smith @*/
6947087cfbeSBarry Smith PetscErrorCode  PetscOptionsInt(const char opt[],const char text[],const char man[],PetscInt defaultv,PetscInt *value,PetscBool  *set)
69553acd3b1SBarry Smith {
69653acd3b1SBarry Smith   PetscErrorCode ierr;
6976356e834SBarry Smith   PetscOptions   amsopt;
69853acd3b1SBarry Smith 
69953acd3b1SBarry Smith   PetscFunctionBegin;
700b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
7016356e834SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_INT,&amsopt);CHKERRQ(ierr);
7026356e834SBarry Smith     ierr = PetscMalloc(sizeof(PetscInt),&amsopt->data);CHKERRQ(ierr);
703a297a907SKarl Rupp 
7046356e834SBarry Smith     *(PetscInt*)amsopt->data = defaultv;
705af6d86caSBarry Smith   }
70653acd3b1SBarry Smith   ierr = PetscOptionsGetInt(PetscOptionsObject.prefix,opt,value,set);CHKERRQ(ierr);
70761b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
7082aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%d>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,defaultv,text,ManSection(man));CHKERRQ(ierr);
70953acd3b1SBarry Smith   }
71053acd3b1SBarry Smith   PetscFunctionReturn(0);
71153acd3b1SBarry Smith }
71253acd3b1SBarry Smith 
71353acd3b1SBarry Smith #undef __FUNCT__
71453acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsString"
71553acd3b1SBarry Smith /*@C
71653acd3b1SBarry Smith    PetscOptionsString - Gets the string value for a particular option in the database.
71753acd3b1SBarry Smith 
7183f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
71953acd3b1SBarry Smith 
72053acd3b1SBarry Smith    Input Parameters:
72153acd3b1SBarry Smith +  opt - option name
72253acd3b1SBarry Smith .  text - short string that describes the option
72353acd3b1SBarry Smith .  man - manual page with additional information on option
724bcbf2dc5SJed Brown .  defaultv - the default (current) value
725bcbf2dc5SJed Brown -  len - length of the result string including null terminator
72653acd3b1SBarry Smith 
72753acd3b1SBarry Smith    Output Parameter:
72853acd3b1SBarry Smith +  value - the value to return
72953acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
73053acd3b1SBarry Smith 
73153acd3b1SBarry Smith    Level: beginner
73253acd3b1SBarry Smith 
73353acd3b1SBarry Smith    Concepts: options database^has int
73453acd3b1SBarry Smith 
73553acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
73653acd3b1SBarry Smith 
7377fccdfe4SBarry 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).
7387fccdfe4SBarry Smith 
73953acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
740acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
741acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsReal(), PetscOptionsBool(),
74253acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
74353acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
744acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
745a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
74653acd3b1SBarry Smith @*/
7477087cfbeSBarry Smith PetscErrorCode  PetscOptionsString(const char opt[],const char text[],const char man[],const char defaultv[],char value[],size_t len,PetscBool  *set)
74853acd3b1SBarry Smith {
74953acd3b1SBarry Smith   PetscErrorCode ierr;
750aee2cecaSBarry Smith   PetscOptions   amsopt;
75153acd3b1SBarry Smith 
75253acd3b1SBarry Smith   PetscFunctionBegin;
753b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
754aee2cecaSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
75564facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
7565b02f95dSBarry Smith     ierr = PetscStrdup(defaultv ? defaultv : "",(char**)&amsopt->data);CHKERRQ(ierr);
757af6d86caSBarry Smith   }
75853acd3b1SBarry Smith   ierr = PetscOptionsGetString(PetscOptionsObject.prefix,opt,value,len,set);CHKERRQ(ierr);
75961b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
7602aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%s>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,defaultv,text,ManSection(man));CHKERRQ(ierr);
76153acd3b1SBarry Smith   }
76253acd3b1SBarry Smith   PetscFunctionReturn(0);
76353acd3b1SBarry Smith }
76453acd3b1SBarry Smith 
76553acd3b1SBarry Smith #undef __FUNCT__
76653acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsReal"
76753acd3b1SBarry Smith /*@C
76853acd3b1SBarry Smith    PetscOptionsReal - Gets the PetscReal value for a particular option in the database.
76953acd3b1SBarry Smith 
7703f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
77153acd3b1SBarry Smith 
77253acd3b1SBarry Smith    Input Parameters:
77353acd3b1SBarry Smith +  opt - option name
77453acd3b1SBarry Smith .  text - short string that describes the option
77553acd3b1SBarry Smith .  man - manual page with additional information on option
77653acd3b1SBarry Smith -  defaultv - the default (current) value
77753acd3b1SBarry Smith 
77853acd3b1SBarry Smith    Output Parameter:
77953acd3b1SBarry Smith +  value - the value to return
78053acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
78153acd3b1SBarry Smith 
78253acd3b1SBarry Smith    Level: beginner
78353acd3b1SBarry Smith 
78453acd3b1SBarry Smith    Concepts: options database^has int
78553acd3b1SBarry Smith 
78653acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
78753acd3b1SBarry Smith 
78853acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
789acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
790acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
79153acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
79253acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
793acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
794a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
79553acd3b1SBarry Smith @*/
7967087cfbeSBarry Smith PetscErrorCode  PetscOptionsReal(const char opt[],const char text[],const char man[],PetscReal defaultv,PetscReal *value,PetscBool  *set)
79753acd3b1SBarry Smith {
79853acd3b1SBarry Smith   PetscErrorCode ierr;
799538aa990SBarry Smith   PetscOptions   amsopt;
80053acd3b1SBarry Smith 
80153acd3b1SBarry Smith   PetscFunctionBegin;
802b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
803538aa990SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_REAL,&amsopt);CHKERRQ(ierr);
804538aa990SBarry Smith     ierr = PetscMalloc(sizeof(PetscReal),&amsopt->data);CHKERRQ(ierr);
805a297a907SKarl Rupp 
806538aa990SBarry Smith     *(PetscReal*)amsopt->data = defaultv;
807538aa990SBarry Smith   }
80853acd3b1SBarry Smith   ierr = PetscOptionsGetReal(PetscOptionsObject.prefix,opt,value,set);CHKERRQ(ierr);
80961b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
8102aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%G>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,defaultv,text,ManSection(man));CHKERRQ(ierr);
81153acd3b1SBarry Smith   }
81253acd3b1SBarry Smith   PetscFunctionReturn(0);
81353acd3b1SBarry Smith }
81453acd3b1SBarry Smith 
81553acd3b1SBarry Smith #undef __FUNCT__
81653acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsScalar"
81753acd3b1SBarry Smith /*@C
81853acd3b1SBarry Smith    PetscOptionsScalar - Gets the scalar value for a particular option in the database.
81953acd3b1SBarry Smith 
8203f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
82153acd3b1SBarry Smith 
82253acd3b1SBarry Smith    Input Parameters:
82353acd3b1SBarry Smith +  opt - option name
82453acd3b1SBarry Smith .  text - short string that describes the option
82553acd3b1SBarry Smith .  man - manual page with additional information on option
82653acd3b1SBarry Smith -  defaultv - the default (current) value
82753acd3b1SBarry Smith 
82853acd3b1SBarry Smith    Output Parameter:
82953acd3b1SBarry Smith +  value - the value to return
83053acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
83153acd3b1SBarry Smith 
83253acd3b1SBarry Smith    Level: beginner
83353acd3b1SBarry Smith 
83453acd3b1SBarry Smith    Concepts: options database^has int
83553acd3b1SBarry Smith 
83653acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
83753acd3b1SBarry Smith 
83853acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
839acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
840acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
84153acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
84253acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
843acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
844a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
84553acd3b1SBarry Smith @*/
8467087cfbeSBarry Smith PetscErrorCode  PetscOptionsScalar(const char opt[],const char text[],const char man[],PetscScalar defaultv,PetscScalar *value,PetscBool  *set)
84753acd3b1SBarry Smith {
84853acd3b1SBarry Smith   PetscErrorCode ierr;
84953acd3b1SBarry Smith 
85053acd3b1SBarry Smith   PetscFunctionBegin;
85153acd3b1SBarry Smith #if !defined(PETSC_USE_COMPLEX)
85253acd3b1SBarry Smith   ierr = PetscOptionsReal(opt,text,man,defaultv,value,set);CHKERRQ(ierr);
85353acd3b1SBarry Smith #else
85453acd3b1SBarry Smith   ierr = PetscOptionsGetScalar(PetscOptionsObject.prefix,opt,value,set);CHKERRQ(ierr);
85553acd3b1SBarry Smith #endif
85653acd3b1SBarry Smith   PetscFunctionReturn(0);
85753acd3b1SBarry Smith }
85853acd3b1SBarry Smith 
85953acd3b1SBarry Smith #undef __FUNCT__
86053acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsName"
86153acd3b1SBarry Smith /*@C
86290d69ab7SBarry 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
86390d69ab7SBarry Smith                       its value is set to false.
86453acd3b1SBarry Smith 
8653f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
86653acd3b1SBarry Smith 
86753acd3b1SBarry Smith    Input Parameters:
86853acd3b1SBarry Smith +  opt - option name
86953acd3b1SBarry Smith .  text - short string that describes the option
87053acd3b1SBarry Smith -  man - manual page with additional information on option
87153acd3b1SBarry Smith 
87253acd3b1SBarry Smith    Output Parameter:
87353acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
87453acd3b1SBarry Smith 
87553acd3b1SBarry Smith    Level: beginner
87653acd3b1SBarry Smith 
87753acd3b1SBarry Smith    Concepts: options database^has int
87853acd3b1SBarry Smith 
87953acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
88053acd3b1SBarry Smith 
88153acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
882acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
883acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
88453acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
88553acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
886acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
887a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
88853acd3b1SBarry Smith @*/
8897087cfbeSBarry Smith PetscErrorCode  PetscOptionsName(const char opt[],const char text[],const char man[],PetscBool  *flg)
89053acd3b1SBarry Smith {
89153acd3b1SBarry Smith   PetscErrorCode ierr;
8921ae3d29cSBarry Smith   PetscOptions   amsopt;
89353acd3b1SBarry Smith 
89453acd3b1SBarry Smith   PetscFunctionBegin;
8951ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
8967781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
897ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
898a297a907SKarl Rupp 
899ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
9001ae3d29cSBarry Smith   }
90153acd3b1SBarry Smith   ierr = PetscOptionsHasName(PetscOptionsObject.prefix,opt,flg);CHKERRQ(ierr);
90261b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
9032aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
90453acd3b1SBarry Smith   }
90553acd3b1SBarry Smith   PetscFunctionReturn(0);
90653acd3b1SBarry Smith }
90753acd3b1SBarry Smith 
90853acd3b1SBarry Smith #undef __FUNCT__
909a264d7a6SBarry Smith #define __FUNCT__ "PetscOptionsFList"
91053acd3b1SBarry Smith /*@C
911a264d7a6SBarry Smith      PetscOptionsFList - Puts a list of option values that a single one may be selected from
91253acd3b1SBarry Smith 
9133f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
91453acd3b1SBarry Smith 
91553acd3b1SBarry Smith    Input Parameters:
91653acd3b1SBarry Smith +  opt - option name
91753acd3b1SBarry Smith .  text - short string that describes the option
91853acd3b1SBarry Smith .  man - manual page with additional information on option
91953acd3b1SBarry Smith .  list - the possible choices
9203cc1e11dSBarry Smith .  defaultv - the default (current) value
9213cc1e11dSBarry Smith -  len - the length of the character array value
92253acd3b1SBarry Smith 
92353acd3b1SBarry Smith    Output Parameter:
92453acd3b1SBarry Smith +  value - the value to return
92553acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
92653acd3b1SBarry Smith 
92753acd3b1SBarry Smith    Level: intermediate
92853acd3b1SBarry Smith 
92953acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
93053acd3b1SBarry Smith 
93153acd3b1SBarry Smith    See PetscOptionsEList() for when the choices are given in a string array
93253acd3b1SBarry Smith 
93353acd3b1SBarry Smith    To get a listing of all currently specified options,
93488c29154SBarry Smith     see PetscOptionsView() or PetscOptionsGetAll()
93553acd3b1SBarry Smith 
936eabe10d7SBarry Smith    Developer Note: This cannot check for invalid selection because of things like MATAIJ that are not included in the list
937eabe10d7SBarry Smith 
93853acd3b1SBarry Smith    Concepts: options database^list
93953acd3b1SBarry Smith 
94053acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
941acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
94253acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
94353acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
944acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
945a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsEnum()
94653acd3b1SBarry Smith @*/
947a264d7a6SBarry Smith PetscErrorCode  PetscOptionsFList(const char opt[],const char ltext[],const char man[],PetscFunctionList list,const char defaultv[],char value[],size_t len,PetscBool  *set)
94853acd3b1SBarry Smith {
94953acd3b1SBarry Smith   PetscErrorCode ierr;
9503cc1e11dSBarry Smith   PetscOptions   amsopt;
95153acd3b1SBarry Smith 
95253acd3b1SBarry Smith   PetscFunctionBegin;
953b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
954a264d7a6SBarry Smith     ierr = PetscOptionsCreate_Private(opt,ltext,man,OPTION_FLIST,&amsopt);CHKERRQ(ierr);
95564facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
9565b02f95dSBarry Smith     ierr = PetscStrdup(defaultv ? defaultv : "",(char**)&amsopt->data);CHKERRQ(ierr);
9573cc1e11dSBarry Smith     amsopt->flist = list;
9583cc1e11dSBarry Smith   }
95953acd3b1SBarry Smith   ierr = PetscOptionsGetString(PetscOptionsObject.prefix,opt,value,len,set);CHKERRQ(ierr);
96061b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
961140e18c1SBarry Smith     ierr = PetscFunctionListPrintTypes(PetscOptionsObject.comm,stdout,PetscOptionsObject.prefix,opt,ltext,man,list,defaultv);CHKERRQ(ierr);CHKERRQ(ierr);
96253acd3b1SBarry Smith   }
96353acd3b1SBarry Smith   PetscFunctionReturn(0);
96453acd3b1SBarry Smith }
96553acd3b1SBarry Smith 
96653acd3b1SBarry Smith #undef __FUNCT__
96753acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEList"
96853acd3b1SBarry Smith /*@C
96953acd3b1SBarry Smith      PetscOptionsEList - Puts a list of option values that a single one may be selected from
97053acd3b1SBarry Smith 
9713f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
97253acd3b1SBarry Smith 
97353acd3b1SBarry Smith    Input Parameters:
97453acd3b1SBarry Smith +  opt - option name
97553acd3b1SBarry Smith .  ltext - short string that describes the option
97653acd3b1SBarry Smith .  man - manual page with additional information on option
977a264d7a6SBarry Smith .  list - the possible choices (one of these must be selected, anything else is invalid)
97853acd3b1SBarry Smith .  ntext - number of choices
97953acd3b1SBarry Smith -  defaultv - the default (current) value
98053acd3b1SBarry Smith 
98153acd3b1SBarry Smith    Output Parameter:
98253acd3b1SBarry Smith +  value - the index of the value to return
98353acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
98453acd3b1SBarry Smith 
98553acd3b1SBarry Smith    Level: intermediate
98653acd3b1SBarry Smith 
98753acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
98853acd3b1SBarry Smith 
989a264d7a6SBarry Smith    See PetscOptionsFList() for when the choices are given in a PetscFunctionList()
99053acd3b1SBarry Smith 
99153acd3b1SBarry Smith    Concepts: options database^list
99253acd3b1SBarry Smith 
99353acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
994acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
99553acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
99653acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
997acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
998a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEnum()
99953acd3b1SBarry Smith @*/
10007087cfbeSBarry 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)
100153acd3b1SBarry Smith {
100253acd3b1SBarry Smith   PetscErrorCode ierr;
100353acd3b1SBarry Smith   PetscInt       i;
10041ae3d29cSBarry Smith   PetscOptions   amsopt;
100553acd3b1SBarry Smith 
100653acd3b1SBarry Smith   PetscFunctionBegin;
10071ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
1008d5649816SBarry Smith     ierr = PetscOptionsCreate_Private(opt,ltext,man,OPTION_ELIST,&amsopt);CHKERRQ(ierr);
100964facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
10105b02f95dSBarry Smith     ierr = PetscStrdup(defaultv ? defaultv : "",(char**)&amsopt->data);CHKERRQ(ierr);
10111ae3d29cSBarry Smith     amsopt->list  = list;
10121ae3d29cSBarry Smith     amsopt->nlist = ntext;
10131ae3d29cSBarry Smith   }
101453acd3b1SBarry Smith   ierr = PetscOptionsGetEList(PetscOptionsObject.prefix,opt,list,ntext,value,set);CHKERRQ(ierr);
101561b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
101653acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%s> (choose one of)",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,defaultv);CHKERRQ(ierr);
101753acd3b1SBarry Smith     for (i=0; i<ntext; i++) {
101853acd3b1SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm," %s",list[i]);CHKERRQ(ierr);
101953acd3b1SBarry Smith     }
10202aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm," (%s)\n",ManSection(man));CHKERRQ(ierr);
102153acd3b1SBarry Smith   }
102253acd3b1SBarry Smith   PetscFunctionReturn(0);
102353acd3b1SBarry Smith }
102453acd3b1SBarry Smith 
102553acd3b1SBarry Smith #undef __FUNCT__
1026acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolGroupBegin"
102753acd3b1SBarry Smith /*@C
1028acfcf0e5SJed Brown      PetscOptionsBoolGroupBegin - First in a series of logical queries on the options database for
1029d5649816SBarry Smith        which at most a single value can be true.
103053acd3b1SBarry Smith 
10313f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
103253acd3b1SBarry Smith 
103353acd3b1SBarry Smith    Input Parameters:
103453acd3b1SBarry Smith +  opt - option name
103553acd3b1SBarry Smith .  text - short string that describes the option
103653acd3b1SBarry Smith -  man - manual page with additional information on option
103753acd3b1SBarry Smith 
103853acd3b1SBarry Smith    Output Parameter:
103953acd3b1SBarry Smith .  flg - whether that option was set or not
104053acd3b1SBarry Smith 
104153acd3b1SBarry Smith    Level: intermediate
104253acd3b1SBarry Smith 
104353acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
104453acd3b1SBarry Smith 
1045acfcf0e5SJed Brown    Must be followed by 0 or more PetscOptionsBoolGroup()s and PetscOptionsBoolGroupEnd()
104653acd3b1SBarry Smith 
104753acd3b1SBarry Smith     Concepts: options database^logical group
104853acd3b1SBarry Smith 
104953acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1050acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
105153acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
105253acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1053acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1054a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
105553acd3b1SBarry Smith @*/
10567087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolGroupBegin(const char opt[],const char text[],const char man[],PetscBool  *flg)
105753acd3b1SBarry Smith {
105853acd3b1SBarry Smith   PetscErrorCode ierr;
10591ae3d29cSBarry Smith   PetscOptions   amsopt;
106053acd3b1SBarry Smith 
106153acd3b1SBarry Smith   PetscFunctionBegin;
10621ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
10637781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1064ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1065a297a907SKarl Rupp 
1066ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
10671ae3d29cSBarry Smith   }
106868b16fdaSBarry Smith   *flg = PETSC_FALSE;
10690298fd71SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,NULL);CHKERRQ(ierr);
107061b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
107153acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  Pick at most one of -------------\n");CHKERRQ(ierr);
10722aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"    -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
107353acd3b1SBarry Smith   }
107453acd3b1SBarry Smith   PetscFunctionReturn(0);
107553acd3b1SBarry Smith }
107653acd3b1SBarry Smith 
107753acd3b1SBarry Smith #undef __FUNCT__
1078acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolGroup"
107953acd3b1SBarry Smith /*@C
1080acfcf0e5SJed Brown      PetscOptionsBoolGroup - One in a series of logical queries on the options database for
1081d5649816SBarry Smith        which at most a single value can be true.
108253acd3b1SBarry Smith 
10833f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
108453acd3b1SBarry Smith 
108553acd3b1SBarry Smith    Input Parameters:
108653acd3b1SBarry Smith +  opt - option name
108753acd3b1SBarry Smith .  text - short string that describes the option
108853acd3b1SBarry Smith -  man - manual page with additional information on option
108953acd3b1SBarry Smith 
109053acd3b1SBarry Smith    Output Parameter:
109153acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
109253acd3b1SBarry Smith 
109353acd3b1SBarry Smith    Level: intermediate
109453acd3b1SBarry Smith 
109553acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
109653acd3b1SBarry Smith 
1097acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin() and preceded a PetscOptionsBoolGroupEnd()
109853acd3b1SBarry Smith 
109953acd3b1SBarry Smith     Concepts: options database^logical group
110053acd3b1SBarry Smith 
110153acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1102acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
110353acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
110453acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1105acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1106a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
110753acd3b1SBarry Smith @*/
11087087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolGroup(const char opt[],const char text[],const char man[],PetscBool  *flg)
110953acd3b1SBarry Smith {
111053acd3b1SBarry Smith   PetscErrorCode ierr;
11111ae3d29cSBarry Smith   PetscOptions   amsopt;
111253acd3b1SBarry Smith 
111353acd3b1SBarry Smith   PetscFunctionBegin;
11141ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
11157781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1116ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1117a297a907SKarl Rupp 
1118ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
11191ae3d29cSBarry Smith   }
112017326d04SJed Brown   *flg = PETSC_FALSE;
11210298fd71SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,NULL);CHKERRQ(ierr);
112261b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
11232aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"    -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
112453acd3b1SBarry Smith   }
112553acd3b1SBarry Smith   PetscFunctionReturn(0);
112653acd3b1SBarry Smith }
112753acd3b1SBarry Smith 
112853acd3b1SBarry Smith #undef __FUNCT__
1129acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolGroupEnd"
113053acd3b1SBarry Smith /*@C
1131acfcf0e5SJed Brown      PetscOptionsBoolGroupEnd - Last in a series of logical queries on the options database for
1132d5649816SBarry Smith        which at most a single value can be true.
113353acd3b1SBarry Smith 
11343f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
113553acd3b1SBarry Smith 
113653acd3b1SBarry Smith    Input Parameters:
113753acd3b1SBarry Smith +  opt - option name
113853acd3b1SBarry Smith .  text - short string that describes the option
113953acd3b1SBarry Smith -  man - manual page with additional information on option
114053acd3b1SBarry Smith 
114153acd3b1SBarry Smith    Output Parameter:
114253acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
114353acd3b1SBarry Smith 
114453acd3b1SBarry Smith    Level: intermediate
114553acd3b1SBarry Smith 
114653acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
114753acd3b1SBarry Smith 
1148acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin()
114953acd3b1SBarry Smith 
115053acd3b1SBarry Smith     Concepts: options database^logical group
115153acd3b1SBarry Smith 
115253acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1153acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
115453acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
115553acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1156acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1157a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
115853acd3b1SBarry Smith @*/
11597087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolGroupEnd(const char opt[],const char text[],const char man[],PetscBool  *flg)
116053acd3b1SBarry Smith {
116153acd3b1SBarry Smith   PetscErrorCode ierr;
11621ae3d29cSBarry Smith   PetscOptions   amsopt;
116353acd3b1SBarry Smith 
116453acd3b1SBarry Smith   PetscFunctionBegin;
11651ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
11667781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1167ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1168a297a907SKarl Rupp 
1169ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
11701ae3d29cSBarry Smith   }
117117326d04SJed Brown   *flg = PETSC_FALSE;
11720298fd71SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,NULL);CHKERRQ(ierr);
117361b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
11742aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"    -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
117553acd3b1SBarry Smith   }
117653acd3b1SBarry Smith   PetscFunctionReturn(0);
117753acd3b1SBarry Smith }
117853acd3b1SBarry Smith 
117953acd3b1SBarry Smith #undef __FUNCT__
1180acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBool"
118153acd3b1SBarry Smith /*@C
1182acfcf0e5SJed Brown    PetscOptionsBool - Determines if a particular option is in the database with a true or false
118353acd3b1SBarry Smith 
11843f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
118553acd3b1SBarry Smith 
118653acd3b1SBarry Smith    Input Parameters:
118753acd3b1SBarry Smith +  opt - option name
118853acd3b1SBarry Smith .  text - short string that describes the option
118953acd3b1SBarry Smith -  man - manual page with additional information on option
119053acd3b1SBarry Smith 
119153acd3b1SBarry Smith    Output Parameter:
119253acd3b1SBarry Smith .  flg - PETSC_TRUE or PETSC_FALSE
119353acd3b1SBarry Smith .  set - PETSC_TRUE if found, else PETSC_FALSE
119453acd3b1SBarry Smith 
119553acd3b1SBarry Smith    Level: beginner
119653acd3b1SBarry Smith 
119753acd3b1SBarry Smith    Concepts: options database^logical
119853acd3b1SBarry Smith 
119953acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
120053acd3b1SBarry Smith 
120153acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
1202acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
1203acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
120453acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
120553acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1206acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1207a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
120853acd3b1SBarry Smith @*/
12097087cfbeSBarry Smith PetscErrorCode  PetscOptionsBool(const char opt[],const char text[],const char man[],PetscBool deflt,PetscBool  *flg,PetscBool  *set)
121053acd3b1SBarry Smith {
121153acd3b1SBarry Smith   PetscErrorCode ierr;
1212ace3abfcSBarry Smith   PetscBool      iset;
1213aee2cecaSBarry Smith   PetscOptions   amsopt;
121453acd3b1SBarry Smith 
121553acd3b1SBarry Smith   PetscFunctionBegin;
1216b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
12177781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1218ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1219a297a907SKarl Rupp 
1220ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = deflt;
1221af6d86caSBarry Smith   }
1222acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,&iset);CHKERRQ(ierr);
122353acd3b1SBarry Smith   if (!iset) {
122453acd3b1SBarry Smith     if (flg) *flg = deflt;
122553acd3b1SBarry Smith   }
122653acd3b1SBarry Smith   if (set) *set = iset;
122761b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
1228ace3abfcSBarry Smith     const char *v = PetscBools[deflt];
12292aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s: <%s> %s (%s)\n",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,v,text,ManSection(man));CHKERRQ(ierr);
123053acd3b1SBarry Smith   }
123153acd3b1SBarry Smith   PetscFunctionReturn(0);
123253acd3b1SBarry Smith }
123353acd3b1SBarry Smith 
123453acd3b1SBarry Smith #undef __FUNCT__
123553acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsRealArray"
123653acd3b1SBarry Smith /*@C
123753acd3b1SBarry Smith    PetscOptionsRealArray - Gets an array of double values for a particular
123853acd3b1SBarry Smith    option in the database. The values must be separated with commas with
123953acd3b1SBarry Smith    no intervening spaces.
124053acd3b1SBarry Smith 
12413f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
124253acd3b1SBarry Smith 
124353acd3b1SBarry Smith    Input Parameters:
124453acd3b1SBarry Smith +  opt - the option one is seeking
124553acd3b1SBarry Smith .  text - short string describing option
124653acd3b1SBarry Smith .  man - manual page for option
124753acd3b1SBarry Smith -  nmax - maximum number of values
124853acd3b1SBarry Smith 
124953acd3b1SBarry Smith    Output Parameter:
125053acd3b1SBarry Smith +  value - location to copy values
125153acd3b1SBarry Smith .  nmax - actual number of values found
125253acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
125353acd3b1SBarry Smith 
125453acd3b1SBarry Smith    Level: beginner
125553acd3b1SBarry Smith 
125653acd3b1SBarry Smith    Notes:
125753acd3b1SBarry Smith    The user should pass in an array of doubles
125853acd3b1SBarry Smith 
125953acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
126053acd3b1SBarry Smith 
126153acd3b1SBarry Smith    Concepts: options database^array of strings
126253acd3b1SBarry Smith 
126353acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1264acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
126553acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
126653acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1267acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1268a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
126953acd3b1SBarry Smith @*/
12707087cfbeSBarry Smith PetscErrorCode  PetscOptionsRealArray(const char opt[],const char text[],const char man[],PetscReal value[],PetscInt *n,PetscBool  *set)
127153acd3b1SBarry Smith {
127253acd3b1SBarry Smith   PetscErrorCode ierr;
127353acd3b1SBarry Smith   PetscInt       i;
1274e26ddf31SBarry Smith   PetscOptions   amsopt;
127553acd3b1SBarry Smith 
127653acd3b1SBarry Smith   PetscFunctionBegin;
1277b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
1278e26ddf31SBarry Smith     PetscReal *vals;
1279e26ddf31SBarry Smith 
1280e26ddf31SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_REAL_ARRAY,&amsopt);CHKERRQ(ierr);
12816e02237eSPeter Brune     ierr = PetscMalloc1((*n),(PetscReal**)&amsopt->data);CHKERRQ(ierr);
1282e26ddf31SBarry Smith     vals = (PetscReal*)amsopt->data;
1283e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1284e26ddf31SBarry Smith     amsopt->arraylength = *n;
1285e26ddf31SBarry Smith   }
128653acd3b1SBarry Smith   ierr = PetscOptionsGetRealArray(PetscOptionsObject.prefix,opt,value,n,set);CHKERRQ(ierr);
128761b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
1288a83599f4SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%G",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,value[0]);CHKERRQ(ierr);
128953acd3b1SBarry Smith     for (i=1; i<*n; i++) {
1290a83599f4SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,",%G",value[i]);CHKERRQ(ierr);
129153acd3b1SBarry Smith     }
12922aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
129353acd3b1SBarry Smith   }
129453acd3b1SBarry Smith   PetscFunctionReturn(0);
129553acd3b1SBarry Smith }
129653acd3b1SBarry Smith 
129753acd3b1SBarry Smith 
129853acd3b1SBarry Smith #undef __FUNCT__
129953acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsIntArray"
130053acd3b1SBarry Smith /*@C
130153acd3b1SBarry Smith    PetscOptionsIntArray - Gets an array of integers for a particular
1302b32a342fSShri Abhyankar    option in the database.
130353acd3b1SBarry Smith 
13043f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
130553acd3b1SBarry Smith 
130653acd3b1SBarry Smith    Input Parameters:
130753acd3b1SBarry Smith +  opt - the option one is seeking
130853acd3b1SBarry Smith .  text - short string describing option
130953acd3b1SBarry Smith .  man - manual page for option
1310f8a50e2bSBarry Smith -  n - maximum number of values
131153acd3b1SBarry Smith 
131253acd3b1SBarry Smith    Output Parameter:
131353acd3b1SBarry Smith +  value - location to copy values
1314f8a50e2bSBarry Smith .  n - actual number of values found
131553acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
131653acd3b1SBarry Smith 
131753acd3b1SBarry Smith    Level: beginner
131853acd3b1SBarry Smith 
131953acd3b1SBarry Smith    Notes:
1320b32a342fSShri Abhyankar    The array can be passed as
1321b32a342fSShri Abhyankar    a comma seperated list:                                 0,1,2,3,4,5,6,7
13220fd488f5SShri Abhyankar    a range (start-end+1):                                  0-8
13230fd488f5SShri Abhyankar    a range with given increment (start-end+1:inc):         0-7:2
13240fd488f5SShri Abhyankar    a combination of values and ranges seperated by commas: 0,1-8,8-15:2
1325b32a342fSShri Abhyankar 
1326b32a342fSShri Abhyankar    There must be no intervening spaces between the values.
132753acd3b1SBarry Smith 
132853acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
132953acd3b1SBarry Smith 
1330b32a342fSShri Abhyankar    Concepts: options database^array of ints
133153acd3b1SBarry Smith 
133253acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1333acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
133453acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
133553acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1336acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1337a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsRealArray()
133853acd3b1SBarry Smith @*/
13397087cfbeSBarry Smith PetscErrorCode  PetscOptionsIntArray(const char opt[],const char text[],const char man[],PetscInt value[],PetscInt *n,PetscBool  *set)
134053acd3b1SBarry Smith {
134153acd3b1SBarry Smith   PetscErrorCode ierr;
134253acd3b1SBarry Smith   PetscInt       i;
1343e26ddf31SBarry Smith   PetscOptions   amsopt;
134453acd3b1SBarry Smith 
134553acd3b1SBarry Smith   PetscFunctionBegin;
1346b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
1347e26ddf31SBarry Smith     PetscInt *vals;
1348e26ddf31SBarry Smith 
1349e26ddf31SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_INT_ARRAY,&amsopt);CHKERRQ(ierr);
13506e02237eSPeter Brune     ierr = PetscMalloc1((*n),(PetscInt**)&amsopt->data);CHKERRQ(ierr);
1351e26ddf31SBarry Smith     vals = (PetscInt*)amsopt->data;
1352e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1353e26ddf31SBarry Smith     amsopt->arraylength = *n;
1354e26ddf31SBarry Smith   }
135553acd3b1SBarry Smith   ierr = PetscOptionsGetIntArray(PetscOptionsObject.prefix,opt,value,n,set);CHKERRQ(ierr);
135661b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
135753acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%d",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,value[0]);CHKERRQ(ierr);
135853acd3b1SBarry Smith     for (i=1; i<*n; i++) {
135953acd3b1SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,",%d",value[i]);CHKERRQ(ierr);
136053acd3b1SBarry Smith     }
13612aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
136253acd3b1SBarry Smith   }
136353acd3b1SBarry Smith   PetscFunctionReturn(0);
136453acd3b1SBarry Smith }
136553acd3b1SBarry Smith 
136653acd3b1SBarry Smith #undef __FUNCT__
136753acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsStringArray"
136853acd3b1SBarry Smith /*@C
136953acd3b1SBarry Smith    PetscOptionsStringArray - Gets an array of string values for a particular
137053acd3b1SBarry Smith    option in the database. The values must be separated with commas with
137153acd3b1SBarry Smith    no intervening spaces.
137253acd3b1SBarry Smith 
13733f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
137453acd3b1SBarry Smith 
137553acd3b1SBarry Smith    Input Parameters:
137653acd3b1SBarry Smith +  opt - the option one is seeking
137753acd3b1SBarry Smith .  text - short string describing option
137853acd3b1SBarry Smith .  man - manual page for option
137953acd3b1SBarry Smith -  nmax - maximum number of strings
138053acd3b1SBarry Smith 
138153acd3b1SBarry Smith    Output Parameter:
138253acd3b1SBarry Smith +  value - location to copy strings
138353acd3b1SBarry Smith .  nmax - actual number of strings found
138453acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
138553acd3b1SBarry Smith 
138653acd3b1SBarry Smith    Level: beginner
138753acd3b1SBarry Smith 
138853acd3b1SBarry Smith    Notes:
138953acd3b1SBarry Smith    The user should pass in an array of pointers to char, to hold all the
139053acd3b1SBarry Smith    strings returned by this function.
139153acd3b1SBarry Smith 
139253acd3b1SBarry Smith    The user is responsible for deallocating the strings that are
139353acd3b1SBarry Smith    returned. The Fortran interface for this routine is not supported.
139453acd3b1SBarry Smith 
139553acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
139653acd3b1SBarry Smith 
139753acd3b1SBarry Smith    Concepts: options database^array of strings
139853acd3b1SBarry Smith 
139953acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1400acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
140153acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
140253acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1403acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1404a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
140553acd3b1SBarry Smith @*/
14067087cfbeSBarry Smith PetscErrorCode  PetscOptionsStringArray(const char opt[],const char text[],const char man[],char *value[],PetscInt *nmax,PetscBool  *set)
140753acd3b1SBarry Smith {
140853acd3b1SBarry Smith   PetscErrorCode ierr;
14091ae3d29cSBarry Smith   PetscOptions   amsopt;
141053acd3b1SBarry Smith 
141153acd3b1SBarry Smith   PetscFunctionBegin;
14121ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
14131ae3d29cSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_STRING_ARRAY,&amsopt);CHKERRQ(ierr);
14146e02237eSPeter Brune     ierr = PetscMalloc1((*nmax),(char**)&amsopt->data);CHKERRQ(ierr);
1415a297a907SKarl Rupp 
14161ae3d29cSBarry Smith     amsopt->arraylength = *nmax;
14171ae3d29cSBarry Smith   }
141853acd3b1SBarry Smith   ierr = PetscOptionsGetStringArray(PetscOptionsObject.prefix,opt,value,nmax,set);CHKERRQ(ierr);
141961b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
14202aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <string1,string2,...>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
142153acd3b1SBarry Smith   }
142253acd3b1SBarry Smith   PetscFunctionReturn(0);
142353acd3b1SBarry Smith }
142453acd3b1SBarry Smith 
1425e2446a98SMatthew Knepley #undef __FUNCT__
1426acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolArray"
1427e2446a98SMatthew Knepley /*@C
1428acfcf0e5SJed Brown    PetscOptionsBoolArray - Gets an array of logical values (true or false) for a particular
1429e2446a98SMatthew Knepley    option in the database. The values must be separated with commas with
1430e2446a98SMatthew Knepley    no intervening spaces.
1431e2446a98SMatthew Knepley 
14323f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
1433e2446a98SMatthew Knepley 
1434e2446a98SMatthew Knepley    Input Parameters:
1435e2446a98SMatthew Knepley +  opt - the option one is seeking
1436e2446a98SMatthew Knepley .  text - short string describing option
1437e2446a98SMatthew Knepley .  man - manual page for option
1438e2446a98SMatthew Knepley -  nmax - maximum number of values
1439e2446a98SMatthew Knepley 
1440e2446a98SMatthew Knepley    Output Parameter:
1441e2446a98SMatthew Knepley +  value - location to copy values
1442e2446a98SMatthew Knepley .  nmax - actual number of values found
1443e2446a98SMatthew Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
1444e2446a98SMatthew Knepley 
1445e2446a98SMatthew Knepley    Level: beginner
1446e2446a98SMatthew Knepley 
1447e2446a98SMatthew Knepley    Notes:
1448e2446a98SMatthew Knepley    The user should pass in an array of doubles
1449e2446a98SMatthew Knepley 
1450e2446a98SMatthew Knepley    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
1451e2446a98SMatthew Knepley 
1452e2446a98SMatthew Knepley    Concepts: options database^array of strings
1453e2446a98SMatthew Knepley 
1454e2446a98SMatthew Knepley .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1455acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
1456e2446a98SMatthew Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1457e2446a98SMatthew Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1458acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1459a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
1460e2446a98SMatthew Knepley @*/
14617087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolArray(const char opt[],const char text[],const char man[],PetscBool value[],PetscInt *n,PetscBool *set)
1462e2446a98SMatthew Knepley {
1463e2446a98SMatthew Knepley   PetscErrorCode ierr;
1464e2446a98SMatthew Knepley   PetscInt       i;
14651ae3d29cSBarry Smith   PetscOptions   amsopt;
1466e2446a98SMatthew Knepley 
1467e2446a98SMatthew Knepley   PetscFunctionBegin;
14681ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
1469ace3abfcSBarry Smith     PetscBool *vals;
14701ae3d29cSBarry Smith 
14717781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL_ARRAY,&amsopt);CHKERRQ(ierr);
14726e02237eSPeter Brune     ierr = PetscMalloc1((*n),(PetscBool**)&amsopt->data);CHKERRQ(ierr);
1473ace3abfcSBarry Smith     vals = (PetscBool*)amsopt->data;
14741ae3d29cSBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
14751ae3d29cSBarry Smith     amsopt->arraylength = *n;
14761ae3d29cSBarry Smith   }
1477acfcf0e5SJed Brown   ierr = PetscOptionsGetBoolArray(PetscOptionsObject.prefix,opt,value,n,set);CHKERRQ(ierr);
1478e2446a98SMatthew Knepley   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
1479e2446a98SMatthew Knepley     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%d",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,value[0]);CHKERRQ(ierr);
1480e2446a98SMatthew Knepley     for (i=1; i<*n; i++) {
1481e2446a98SMatthew Knepley       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,",%d",value[i]);CHKERRQ(ierr);
1482e2446a98SMatthew Knepley     }
14832aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
1484e2446a98SMatthew Knepley   }
1485e2446a98SMatthew Knepley   PetscFunctionReturn(0);
1486e2446a98SMatthew Knepley }
1487e2446a98SMatthew Knepley 
14888cc676e6SMatthew G Knepley #undef __FUNCT__
14898cc676e6SMatthew G Knepley #define __FUNCT__ "PetscOptionsViewer"
14908cc676e6SMatthew G Knepley /*@C
1491*d1da0b69SBarry Smith    PetscOptionsViewer - Gets a viewer appropriate for the type indicated by the user
14928cc676e6SMatthew G Knepley 
14938cc676e6SMatthew G Knepley    Logically Collective on the communicator passed in PetscOptionsBegin()
14948cc676e6SMatthew G Knepley 
14958cc676e6SMatthew G Knepley    Input Parameters:
14968cc676e6SMatthew G Knepley +  opt - option name
14978cc676e6SMatthew G Knepley .  text - short string that describes the option
14988cc676e6SMatthew G Knepley -  man - manual page with additional information on option
14998cc676e6SMatthew G Knepley 
15008cc676e6SMatthew G Knepley    Output Parameter:
15018cc676e6SMatthew G Knepley +  viewer - the viewer
15028cc676e6SMatthew G Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
15038cc676e6SMatthew G Knepley 
15048cc676e6SMatthew G Knepley    Level: beginner
15058cc676e6SMatthew G Knepley 
15068cc676e6SMatthew G Knepley    Concepts: options database^has int
15078cc676e6SMatthew G Knepley 
15088cc676e6SMatthew G Knepley    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
15098cc676e6SMatthew G Knepley 
1510*d1da0b69SBarry Smith    See PetscOptionsGetVieweer() for the format of the supplied viewer and its options
15118cc676e6SMatthew G Knepley 
15128cc676e6SMatthew G Knepley .seealso: PetscOptionsGetViewer(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
15138cc676e6SMatthew G Knepley           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
15148cc676e6SMatthew G Knepley           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
15158cc676e6SMatthew G Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
15168cc676e6SMatthew G Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
15178cc676e6SMatthew G Knepley           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1518a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
15198cc676e6SMatthew G Knepley @*/
1520cffb1e40SBarry Smith PetscErrorCode  PetscOptionsViewer(const char opt[],const char text[],const char man[],PetscViewer *viewer,PetscViewerFormat *format,PetscBool  *set)
15218cc676e6SMatthew G Knepley {
15228cc676e6SMatthew G Knepley   PetscErrorCode ierr;
15238cc676e6SMatthew G Knepley   PetscOptions   amsopt;
15248cc676e6SMatthew G Knepley 
15258cc676e6SMatthew G Knepley   PetscFunctionBegin;
15268cc676e6SMatthew G Knepley   if (!PetscOptionsPublishCount) {
15278cc676e6SMatthew G Knepley     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
152864facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
15295b02f95dSBarry Smith     ierr = PetscStrdup("",(char**)&amsopt->data);CHKERRQ(ierr);
15308cc676e6SMatthew G Knepley   }
1531cffb1e40SBarry Smith   ierr = PetscOptionsGetViewer(PetscOptionsObject.comm,PetscOptionsObject.prefix,opt,viewer,format,set);CHKERRQ(ierr);
15328cc676e6SMatthew G Knepley   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
15338cc676e6SMatthew G Knepley     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%s>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,"",text,ManSection(man));CHKERRQ(ierr);
15348cc676e6SMatthew G Knepley   }
15358cc676e6SMatthew G Knepley   PetscFunctionReturn(0);
15368cc676e6SMatthew G Knepley }
15378cc676e6SMatthew G Knepley 
153853acd3b1SBarry Smith 
153953acd3b1SBarry Smith #undef __FUNCT__
154053acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsHead"
154153acd3b1SBarry Smith /*@C
1542b52f573bSBarry Smith      PetscOptionsHead - Puts a heading before listing any more published options. Used, for example,
154353acd3b1SBarry Smith             in KSPSetFromOptions_GMRES().
154453acd3b1SBarry Smith 
15453f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
154653acd3b1SBarry Smith 
154753acd3b1SBarry Smith    Input Parameter:
154853acd3b1SBarry Smith .   head - the heading text
154953acd3b1SBarry Smith 
155053acd3b1SBarry Smith 
155153acd3b1SBarry Smith    Level: intermediate
155253acd3b1SBarry Smith 
155353acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
155453acd3b1SBarry Smith 
1555b52f573bSBarry Smith           Can be followed by a call to PetscOptionsTail() in the same function.
155653acd3b1SBarry Smith 
155753acd3b1SBarry Smith    Concepts: options database^subheading
155853acd3b1SBarry Smith 
155953acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1560acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
156153acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
156253acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1563acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1564a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
156553acd3b1SBarry Smith @*/
15667087cfbeSBarry Smith PetscErrorCode  PetscOptionsHead(const char head[])
156753acd3b1SBarry Smith {
156853acd3b1SBarry Smith   PetscErrorCode ierr;
156953acd3b1SBarry Smith 
157053acd3b1SBarry Smith   PetscFunctionBegin;
157161b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
157253acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  %s\n",head);CHKERRQ(ierr);
157353acd3b1SBarry Smith   }
157453acd3b1SBarry Smith   PetscFunctionReturn(0);
157553acd3b1SBarry Smith }
157653acd3b1SBarry Smith 
157753acd3b1SBarry Smith 
157853acd3b1SBarry Smith 
157953acd3b1SBarry Smith 
158053acd3b1SBarry Smith 
158153acd3b1SBarry Smith 
1582