xref: /petsc/src/sys/objects/aoptions.c (revision cffb1e400d6c3ea2f6cb522ae2432dc42cf29e73)
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 
8c6db04a5SJed Brown #include <petscsys.h>        /*I  "petscsys.h"   I*/
953acd3b1SBarry Smith #if defined(PETSC_HAVE_STDLIB_H)
1053acd3b1SBarry Smith #include <stdlib.h>
1153acd3b1SBarry Smith #endif
1253acd3b1SBarry Smith 
132aa6d131SJed Brown #define ManSection(str) ((str)?(str):"None")
142aa6d131SJed Brown 
1553acd3b1SBarry Smith /*
1653acd3b1SBarry Smith     Keep a linked list of options that have been posted and we are waiting for
173fc1eb6aSBarry Smith    user selection. See the manual page for PetscOptionsBegin()
1853acd3b1SBarry Smith 
1953acd3b1SBarry Smith     Eventually we'll attach this beast to a MPI_Comm
2053acd3b1SBarry Smith */
21f8d0b74dSMatthew Knepley PetscOptionsObjectType PetscOptionsObject;
2253acd3b1SBarry Smith PetscInt               PetscOptionsPublishCount = 0;
2353acd3b1SBarry Smith 
2453acd3b1SBarry Smith #undef __FUNCT__
2553acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsBegin_Private"
2653acd3b1SBarry Smith /*
2753acd3b1SBarry Smith     Handles setting up the data structure in a call to PetscOptionsBegin()
2853acd3b1SBarry Smith */
2953acd3b1SBarry Smith PetscErrorCode PetscOptionsBegin_Private(MPI_Comm comm,const char prefix[],const char title[],const char mansec[])
3053acd3b1SBarry Smith {
3153acd3b1SBarry Smith   PetscErrorCode ierr;
3253acd3b1SBarry Smith 
3353acd3b1SBarry Smith   PetscFunctionBegin;
3453acd3b1SBarry Smith   PetscOptionsObject.next          = 0;
3553acd3b1SBarry Smith   PetscOptionsObject.comm          = comm;
366356e834SBarry Smith   PetscOptionsObject.changedmethod = PETSC_FALSE;
37c31cb41cSBarry Smith   ierr = PetscFree(PetscOptionsObject.prefix);CHKERRQ(ierr);
3853acd3b1SBarry Smith   ierr = PetscStrallocpy(prefix,&PetscOptionsObject.prefix);CHKERRQ(ierr);
39c31cb41cSBarry Smith   ierr = PetscFree(PetscOptionsObject.title);CHKERRQ(ierr);
4053acd3b1SBarry Smith   ierr = PetscStrallocpy(title,&PetscOptionsObject.title);CHKERRQ(ierr);
4153acd3b1SBarry Smith 
4253acd3b1SBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-help",&PetscOptionsObject.printhelp);CHKERRQ(ierr);
4353acd3b1SBarry Smith   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1) {
4461b37b28SSatish Balay     if (!PetscOptionsObject.alreadyprinted) {
4553acd3b1SBarry Smith       ierr = (*PetscHelpPrintf)(comm,"%s -------------------------------------------------\n",title);CHKERRQ(ierr);
4653acd3b1SBarry Smith     }
4761b37b28SSatish Balay   }
4853acd3b1SBarry Smith   PetscFunctionReturn(0);
4953acd3b1SBarry Smith }
5053acd3b1SBarry Smith 
513194b578SJed Brown #undef __FUNCT__
523194b578SJed Brown #define __FUNCT__ "PetscObjectOptionsBegin_Private"
533194b578SJed Brown /*
543194b578SJed Brown     Handles setting up the data structure in a call to PetscObjectOptionsBegin()
553194b578SJed Brown */
563194b578SJed Brown PetscErrorCode PetscObjectOptionsBegin_Private(PetscObject obj)
573194b578SJed Brown {
583194b578SJed Brown   PetscErrorCode ierr;
593194b578SJed Brown   char title[256];
603194b578SJed Brown   PetscBool flg;
613194b578SJed Brown 
623194b578SJed Brown   PetscFunctionBegin;
633194b578SJed Brown   PetscValidHeader(obj,1);
643194b578SJed Brown   PetscOptionsObject.object = obj;
653194b578SJed Brown   PetscOptionsObject.alreadyprinted = obj->optionsprinted;
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 
916e655a9bSJed Brown   ierr             = PetscNew(struct _n_PetscOptions,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 
10153acd3b1SBarry Smith   if (!PetscOptionsObject.next) {
10253acd3b1SBarry Smith     PetscOptionsObject.next = *amsopt;
10353acd3b1SBarry Smith   } else {
10453acd3b1SBarry Smith     next = PetscOptionsObject.next;
10553acd3b1SBarry Smith     while (next->next) next = next->next;
10653acd3b1SBarry Smith     next->next = *amsopt;
10753acd3b1SBarry Smith   }
10853acd3b1SBarry Smith   PetscFunctionReturn(0);
10953acd3b1SBarry Smith }
11053acd3b1SBarry Smith 
11153acd3b1SBarry Smith #undef __FUNCT__
112aee2cecaSBarry Smith #define __FUNCT__ "PetscScanString"
113aee2cecaSBarry Smith /*
1143fc1eb6aSBarry Smith     PetscScanString -  Gets user input via stdin from process and broadcasts to all processes
1153fc1eb6aSBarry Smith 
1163fc1eb6aSBarry Smith     Collective on MPI_Comm
1173fc1eb6aSBarry Smith 
1183fc1eb6aSBarry Smith    Input Parameters:
1193fc1eb6aSBarry Smith +     commm - communicator for the broadcast, must be PETSC_COMM_WORLD
1203fc1eb6aSBarry Smith .     n - length of the string, must be the same on all processes
1213fc1eb6aSBarry Smith -     str - location to store input
122aee2cecaSBarry Smith 
123aee2cecaSBarry Smith     Bugs:
124aee2cecaSBarry Smith .   Assumes process 0 of the given communicator has access to stdin
125aee2cecaSBarry Smith 
126aee2cecaSBarry Smith */
1273fc1eb6aSBarry Smith static PetscErrorCode PetscScanString(MPI_Comm comm,size_t n,char str[])
128aee2cecaSBarry Smith {
129330cf3c9SBarry Smith   size_t         i;
130aee2cecaSBarry Smith   char           c;
1313fc1eb6aSBarry Smith   PetscMPIInt    rank,nm;
132aee2cecaSBarry Smith   PetscErrorCode ierr;
133aee2cecaSBarry Smith 
134aee2cecaSBarry Smith   PetscFunctionBegin;
135aee2cecaSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
136aee2cecaSBarry Smith   if (!rank) {
137aee2cecaSBarry Smith     c = (char) getchar();
138aee2cecaSBarry Smith     i = 0;
139aee2cecaSBarry Smith     while ( c != '\n' && i < n-1) {
140aee2cecaSBarry Smith       str[i++] = c;
141aee2cecaSBarry Smith       c = (char) getchar();
142aee2cecaSBarry Smith     }
143aee2cecaSBarry Smith     str[i] = 0;
144aee2cecaSBarry Smith   }
1453fc1eb6aSBarry Smith   nm   = PetscMPIIntCast(n);
1463fc1eb6aSBarry Smith   ierr = MPI_Bcast(str,nm,MPI_CHAR,0,comm);CHKERRQ(ierr);
147aee2cecaSBarry Smith   PetscFunctionReturn(0);
148aee2cecaSBarry Smith }
149aee2cecaSBarry Smith 
150aee2cecaSBarry Smith #undef __FUNCT__
151aee2cecaSBarry Smith #define __FUNCT__ "PetscOptionsGetFromTextInput"
152aee2cecaSBarry Smith /*
1533cc1e11dSBarry Smith     PetscOptionsGetFromTextInput - Presents all the PETSc Options processed by the program so the user may change them at runtime
154aee2cecaSBarry Smith 
155aee2cecaSBarry Smith     Notes: this isn't really practical, it is just to demonstrate the principle
156aee2cecaSBarry Smith 
157aee2cecaSBarry Smith     Bugs:
158aee2cecaSBarry Smith +    All processes must traverse through the exact same set of option queries do to the call to PetscScanString()
1593cc1e11dSBarry Smith .    Internal strings have arbitrary length and string copies are not checked that they fit into string space
160aee2cecaSBarry Smith -    Only works for PetscInt == int, PetscReal == double etc
161aee2cecaSBarry Smith 
1623cc1e11dSBarry Smith     Developer Notes: Normally the GUI that presents the options the user and retrieves the values would be running in a different
1633cc1e11dSBarry Smith      address space and communicating with the PETSc program
1643cc1e11dSBarry Smith 
165aee2cecaSBarry Smith */
166aee2cecaSBarry Smith PetscErrorCode PetscOptionsGetFromTextInput()
1676356e834SBarry Smith {
1686356e834SBarry Smith   PetscErrorCode ierr;
1696356e834SBarry Smith   PetscOptions   next = PetscOptionsObject.next;
1706356e834SBarry Smith   char           str[512];
171a4404d99SBarry Smith   PetscInt       id;
172a4404d99SBarry Smith   PetscReal      ir,*valr;
173330cf3c9SBarry Smith   PetscInt       *vald;
174330cf3c9SBarry Smith   size_t         i;
1756356e834SBarry Smith 
176e26ddf31SBarry Smith   ierr = (*PetscPrintf)(PETSC_COMM_WORLD,"%s -------------------------------------------------\n",PetscOptionsObject.title);CHKERRQ(ierr);
1776356e834SBarry Smith   while (next) {
1786356e834SBarry Smith     switch (next->type) {
1796356e834SBarry Smith       case OPTION_HEAD:
1806356e834SBarry Smith         break;
181e26ddf31SBarry Smith       case OPTION_INT_ARRAY:
182e26ddf31SBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",next->option+1);CHKERRQ(ierr);
183e26ddf31SBarry Smith         vald = (PetscInt*) next->data;
184e26ddf31SBarry Smith         for (i=0; i<next->arraylength; i++) {
185e26ddf31SBarry Smith           ierr = PetscPrintf(PETSC_COMM_WORLD,"%d",vald[i]);CHKERRQ(ierr);
186e26ddf31SBarry Smith           if (i < next->arraylength-1) {
187e26ddf31SBarry Smith             ierr = PetscPrintf(PETSC_COMM_WORLD,",");CHKERRQ(ierr);
188e26ddf31SBarry Smith           }
189e26ddf31SBarry Smith         }
190e26ddf31SBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,">: %s (%s)",next->text,next->man);CHKERRQ(ierr);
191e26ddf31SBarry Smith         ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
192e26ddf31SBarry Smith         if (str[0]) {
193e26ddf31SBarry Smith           PetscToken token;
194e26ddf31SBarry Smith           PetscInt   n=0,nmax = next->arraylength,*dvalue = (PetscInt*)next->data,start,end;
195e26ddf31SBarry Smith           size_t     len;
196e26ddf31SBarry Smith           char*      value;
197ace3abfcSBarry Smith           PetscBool  foundrange;
198e26ddf31SBarry Smith 
199e26ddf31SBarry Smith           next->set = PETSC_TRUE;
200e26ddf31SBarry Smith           value = str;
201e26ddf31SBarry Smith 	  ierr = PetscTokenCreate(value,',',&token);CHKERRQ(ierr);
202e26ddf31SBarry Smith 	  ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
203e26ddf31SBarry Smith 	  while (n < nmax) {
204e26ddf31SBarry Smith 	    if (!value) break;
205e26ddf31SBarry Smith 
206e26ddf31SBarry Smith 	    /* look for form  d-D where d and D are integers */
207e26ddf31SBarry Smith 	    foundrange = PETSC_FALSE;
208e26ddf31SBarry Smith 	    ierr      = PetscStrlen(value,&len);CHKERRQ(ierr);
209e26ddf31SBarry Smith 	    if (value[0] == '-') i=2;
210e26ddf31SBarry Smith 	    else i=1;
211330cf3c9SBarry Smith 	    for (;i<len; i++) {
212e26ddf31SBarry Smith 	      if (value[i] == '-') {
213e32f2f54SBarry Smith 		if (i == len-1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Error in %D-th array entry %s\n",n,value);
214e26ddf31SBarry Smith 		value[i] = 0;
215cfbddea1SSatish Balay 		ierr     = PetscOptionsStringToInt(value,&start);CHKERRQ(ierr);
216cfbddea1SSatish Balay 		ierr     = PetscOptionsStringToInt(value+i+1,&end);CHKERRQ(ierr);
217e32f2f54SBarry 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);
218e32f2f54SBarry 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);
219e26ddf31SBarry Smith 		for (;start<end; start++) {
220e26ddf31SBarry Smith 		  *dvalue = start; dvalue++;n++;
221e26ddf31SBarry Smith 		}
222e26ddf31SBarry Smith 		foundrange = PETSC_TRUE;
223e26ddf31SBarry Smith 		break;
224e26ddf31SBarry Smith 	      }
225e26ddf31SBarry Smith 	    }
226e26ddf31SBarry Smith 	    if (!foundrange) {
227cfbddea1SSatish Balay 	      ierr      = PetscOptionsStringToInt(value,dvalue);CHKERRQ(ierr);
228e26ddf31SBarry Smith 	      dvalue++;
229e26ddf31SBarry Smith 	      n++;
230e26ddf31SBarry Smith 	    }
231e26ddf31SBarry Smith 	    ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
232e26ddf31SBarry Smith 	  }
2338c74ee41SBarry Smith 	  ierr = PetscTokenDestroy(&token);CHKERRQ(ierr);
234e26ddf31SBarry Smith         }
235e26ddf31SBarry Smith         break;
236e26ddf31SBarry Smith       case OPTION_REAL_ARRAY:
237e26ddf31SBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",next->option+1);CHKERRQ(ierr);
238e26ddf31SBarry Smith         valr = (PetscReal*) next->data;
239e26ddf31SBarry Smith         for (i=0; i<next->arraylength; i++) {
240e26ddf31SBarry Smith           ierr = PetscPrintf(PETSC_COMM_WORLD,"%g",valr[i]);CHKERRQ(ierr);
241e26ddf31SBarry Smith           if (i < next->arraylength-1) {
242e26ddf31SBarry Smith             ierr = PetscPrintf(PETSC_COMM_WORLD,",");CHKERRQ(ierr);
243e26ddf31SBarry Smith           }
244e26ddf31SBarry Smith         }
245e26ddf31SBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,">: %s (%s)",next->text,next->man);CHKERRQ(ierr);
246e26ddf31SBarry Smith         ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
247e26ddf31SBarry Smith         if (str[0]) {
248e26ddf31SBarry Smith           PetscToken token;
249e26ddf31SBarry Smith           PetscInt   n=0,nmax = next->arraylength;
250e26ddf31SBarry Smith           PetscReal   *dvalue = (PetscReal*)next->data;
251e26ddf31SBarry Smith           char*      value;
252e26ddf31SBarry Smith 
253e26ddf31SBarry Smith           next->set = PETSC_TRUE;
254e26ddf31SBarry Smith           value = str;
255e26ddf31SBarry Smith 	  ierr = PetscTokenCreate(value,',',&token);CHKERRQ(ierr);
256e26ddf31SBarry Smith 	  ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
257e26ddf31SBarry Smith 	  while (n < nmax) {
258e26ddf31SBarry Smith 	    if (!value) break;
259cfbddea1SSatish Balay             ierr      = PetscOptionsStringToReal(value,dvalue);CHKERRQ(ierr);
260e26ddf31SBarry Smith 	    dvalue++;
261e26ddf31SBarry Smith 	    n++;
262e26ddf31SBarry Smith 	    ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
263e26ddf31SBarry Smith 	  }
2648c74ee41SBarry Smith 	  ierr = PetscTokenDestroy(&token);CHKERRQ(ierr);
265e26ddf31SBarry Smith         }
266e26ddf31SBarry Smith         break;
2676356e834SBarry Smith       case OPTION_INT:
268e26ddf31SBarry 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);
2693fc1eb6aSBarry Smith         ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
2703fc1eb6aSBarry Smith         if (str[0]) {
271c272547aSJed Brown #if defined(PETSC_USE_64BIT_INDICES)
272c272547aSJed Brown           sscanf(str,"%lld",&id);
273c272547aSJed Brown #else
274aee2cecaSBarry Smith           sscanf(str,"%d",&id);
275c272547aSJed Brown #endif
276aee2cecaSBarry Smith           next->set = PETSC_TRUE;
277aee2cecaSBarry Smith           *((PetscInt*)next->data) = id;
278aee2cecaSBarry Smith         }
279aee2cecaSBarry Smith         break;
280aee2cecaSBarry Smith       case OPTION_REAL:
281e26ddf31SBarry 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);
2823fc1eb6aSBarry Smith         ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
2833fc1eb6aSBarry Smith         if (str[0]) {
284ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
285a4404d99SBarry Smith           sscanf(str,"%e",&ir);
286ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL_DOUBLE)
287aee2cecaSBarry Smith           sscanf(str,"%le",&ir);
288ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL___FLOAT128)
289d9822059SBarry Smith           ir = strtoflt128(str,0);
290d9822059SBarry Smith #else
291513dbe71SLisandro Dalcin           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Unknown scalar type");
292a4404d99SBarry Smith #endif
293aee2cecaSBarry Smith           next->set = PETSC_TRUE;
294aee2cecaSBarry Smith           *((PetscReal*)next->data) = ir;
295aee2cecaSBarry Smith         }
296aee2cecaSBarry Smith         break;
297aee2cecaSBarry Smith       case OPTION_LOGICAL:
298aee2cecaSBarry Smith       case OPTION_STRING:
299e26ddf31SBarry 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);
3003fc1eb6aSBarry Smith         ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3013fc1eb6aSBarry Smith         if (str[0]) {
302aee2cecaSBarry Smith           next->set = PETSC_TRUE;
303330cf3c9SBarry Smith           ierr = PetscStrcpy((char*)next->data,str);CHKERRQ(ierr);
3046356e834SBarry Smith         }
3056356e834SBarry Smith         break;
3063cc1e11dSBarry Smith       case OPTION_LIST:
3073cc1e11dSBarry Smith         ierr = PetscFListPrintTypes(PETSC_COMM_WORLD,stdout,PetscOptionsObject.prefix,next->option,next->text,next->man,next->flist,(char*)next->data);CHKERRQ(ierr);
3083cc1e11dSBarry Smith         ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3093cc1e11dSBarry Smith         if (str[0]) {
3103cc1e11dSBarry Smith 	  PetscOptionsObject.changedmethod = PETSC_TRUE;
3113cc1e11dSBarry Smith           next->set = PETSC_TRUE;
312330cf3c9SBarry Smith           ierr = PetscStrcpy((char*)next->data,str);CHKERRQ(ierr);
3133cc1e11dSBarry Smith         }
3143cc1e11dSBarry Smith         break;
315b432afdaSMatthew Knepley     default:
316b432afdaSMatthew Knepley       break;
3176356e834SBarry Smith     }
3186356e834SBarry Smith     next = next->next;
3196356e834SBarry Smith   }
3206356e834SBarry Smith   PetscFunctionReturn(0);
3216356e834SBarry Smith }
3226356e834SBarry Smith 
323b3506946SBarry Smith #if defined(PETSC_HAVE_AMS)
3241ae3d29cSBarry Smith #define CHKERRAMS(err)  if (err) {char *msg; AMS_Explain_error((err), &(msg)); SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"AMS Error: %s",msg);}
3251ae3d29cSBarry Smith #define CHKERRAMSFieldName(err,fn)  if (err) {char *msg; AMS_Explain_error((err), &(msg)); SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB,"Fieldname %s, AMS Error: %s",fn,msg);}
326d5649816SBarry Smith 
327d5649816SBarry Smith static int  count = 0;
328d5649816SBarry Smith 
329b3506946SBarry Smith #undef __FUNCT__
330d5649816SBarry Smith #define __FUNCT__ "PetscOptionsAMSDestroy"
331d5649816SBarry Smith PetscErrorCode PetscOptionsAMSDestroy(void)
332d5649816SBarry Smith {
333d5649816SBarry Smith   PetscErrorCode ierr;
334d5649816SBarry Smith   AMS_Comm       acomm = -1;
335d5649816SBarry Smith   AMS_Memory     amem = -1;
336d5649816SBarry Smith   char           options[16];
337d5649816SBarry Smith   const char     *string = "Exit";
338d5649816SBarry Smith 
339d5649816SBarry Smith   /* the next line is a bug, this will only work if all processors are here, the comm passed in is ignored!!! */
340d5649816SBarry Smith   ierr = PetscViewerAMSGetAMSComm(PETSC_VIEWER_AMS_(PETSC_COMM_WORLD),&acomm);CHKERRQ(ierr);
341d5649816SBarry Smith   sprintf(options,"Options_%d",count++);
342d5649816SBarry Smith   ierr = AMS_Memory_create(acomm,options,&amem);CHKERRAMS(ierr);
343d5649816SBarry Smith   ierr = AMS_Memory_add_field(amem,"Exit",&string,1,AMS_STRING,AMS_READ,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,"Exit");
344d5649816SBarry Smith 
345d5649816SBarry Smith   ierr = AMS_Memory_take_access(amem);CHKERRAMS(ierr);
346d5649816SBarry Smith   ierr = AMS_Memory_publish(amem);CHKERRAMS(ierr);
347d5649816SBarry Smith   ierr = AMS_Memory_grant_access(amem);CHKERRAMS(ierr);
348d5649816SBarry Smith   /* wait until accessor has unlocked the memory */
349d5649816SBarry Smith   ierr = AMS_Memory_lock(amem,0);CHKERRAMS(ierr);
350d5649816SBarry Smith   ierr = AMS_Memory_take_access(amem);CHKERRAMS(ierr);
351d5649816SBarry Smith   ierr = AMS_Memory_grant_access(amem);CHKERRAMS(ierr);
352d5649816SBarry Smith   ierr = AMS_Memory_destroy(amem);CHKERRAMS(ierr);
353d5649816SBarry Smith   PetscFunctionReturn(0);
354d5649816SBarry Smith }
355d5649816SBarry Smith 
356d5649816SBarry Smith #undef __FUNCT__
357d5649816SBarry Smith #define __FUNCT__ "PetscOptionsAMSInput"
358b3506946SBarry Smith /*
359d5649816SBarry Smith     PetscOptionsAMSInput - Presents all the PETSc Options processed by the program so the user may change them at runtime using the AMS
360b3506946SBarry Smith 
361b3506946SBarry Smith     Bugs:
362b3506946SBarry Smith +    All processes must traverse through the exact same set of option queries do to the call to PetscScanString()
363b3506946SBarry Smith .    Internal strings have arbitrary length and string copies are not checked that they fit into string space
364b3506946SBarry Smith -    Only works for PetscInt == int, PetscReal == double etc
365b3506946SBarry Smith 
366b3506946SBarry Smith 
367b3506946SBarry Smith */
368d5649816SBarry Smith PetscErrorCode PetscOptionsAMSInput()
369b3506946SBarry Smith {
370b3506946SBarry Smith   PetscErrorCode ierr;
371b3506946SBarry Smith   PetscOptions   next = PetscOptionsObject.next;
372d5649816SBarry Smith   static int     mancount = 0;
373b3506946SBarry Smith   char           options[16];
374b3506946SBarry Smith   AMS_Comm       acomm = -1;
375b3506946SBarry Smith   AMS_Memory     amem = -1;
376ace3abfcSBarry Smith   PetscBool      changedmethod = PETSC_FALSE;
3779f32e415SBarry Smith   char           manname[16];
378b3506946SBarry Smith 
379b3506946SBarry Smith   /* the next line is a bug, this will only work if all processors are here, the comm passed in is ignored!!! */
380b3506946SBarry Smith   ierr = PetscViewerAMSGetAMSComm(PETSC_VIEWER_AMS_(PETSC_COMM_WORLD),&acomm);CHKERRQ(ierr);
381b3506946SBarry Smith   sprintf(options,"Options_%d",count++);
3821ae3d29cSBarry Smith   ierr = AMS_Memory_create(acomm,options,&amem);CHKERRAMS(ierr);
3831ae3d29cSBarry Smith   ierr = AMS_Memory_take_access(amem);CHKERRAMS(ierr);
3841bc75a8dSBarry Smith   PetscOptionsObject.pprefix = PetscOptionsObject.prefix; /* AMS will change this, so cannot pass prefix directly */
3851bc75a8dSBarry Smith 
3861ae3d29cSBarry Smith   ierr = AMS_Memory_add_field(amem,PetscOptionsObject.title,&PetscOptionsObject.pprefix,1,AMS_STRING,AMS_READ,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,PetscOptionsObject.title);
3870c74a584SJed Brown   /* ierr = AMS_Memory_add_field(amem,"mansec",&PetscOptionsObject.pprefix,1,AMS_STRING,AMS_READ,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMS(ierr); */
3881ae3d29cSBarry Smith   ierr = AMS_Memory_add_field(amem,"ChangedMethod",&changedmethod,1,AMS_BOOLEAN,AMS_WRITE,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,"ChangedMethod");
389b3506946SBarry Smith 
390b3506946SBarry Smith   while (next) {
3911ae3d29cSBarry Smith     ierr = AMS_Memory_add_field(amem,next->option,&next->set,1,AMS_INT,AMS_WRITE,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,next->option);
3921bc75a8dSBarry Smith     ierr =  PetscMalloc(sizeof(char*),&next->pman);CHKERRQ(ierr);
3931bc75a8dSBarry Smith     *(char **)next->pman = next->man;
3949f32e415SBarry Smith     sprintf(manname,"man_%d",mancount++);
3951ae3d29cSBarry Smith     ierr = AMS_Memory_add_field(amem,manname,next->pman,1,AMS_STRING,AMS_READ,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,manname);
3969f32e415SBarry Smith 
397b3506946SBarry Smith     switch (next->type) {
398b3506946SBarry Smith       case OPTION_HEAD:
399b3506946SBarry Smith         break;
400b3506946SBarry Smith       case OPTION_INT_ARRAY:
4011ae3d29cSBarry Smith         ierr = AMS_Memory_add_field(amem,next->text,next->data,next->arraylength,AMS_INT,AMS_WRITE,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,next->text);
402b3506946SBarry Smith         break;
403b3506946SBarry Smith       case OPTION_REAL_ARRAY:
4041ae3d29cSBarry Smith         ierr = AMS_Memory_add_field(amem,next->text,next->data,next->arraylength,AMS_DOUBLE,AMS_WRITE,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,next->text);
405b3506946SBarry Smith         break;
406b3506946SBarry Smith       case OPTION_INT:
4071ae3d29cSBarry Smith 	ierr = AMS_Memory_add_field(amem,next->text,next->data,1,AMS_INT,AMS_WRITE,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,next->text);
408b3506946SBarry Smith         break;
409b3506946SBarry Smith       case OPTION_REAL:
4101ae3d29cSBarry Smith 	ierr = AMS_Memory_add_field(amem,next->text,next->data,1,AMS_DOUBLE,AMS_WRITE,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,next->text);
411b3506946SBarry Smith         break;
412b3506946SBarry Smith       case OPTION_LOGICAL:
4131ae3d29cSBarry Smith 	ierr = AMS_Memory_add_field(amem,next->text,next->data,1,AMS_BOOLEAN,AMS_WRITE,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,next->text);
4141ae3d29cSBarry Smith         break;
4151ae3d29cSBarry Smith       case OPTION_LOGICAL_ARRAY:
4161ae3d29cSBarry Smith 	ierr = AMS_Memory_add_field(amem,next->text,next->data,next->arraylength,AMS_BOOLEAN,AMS_WRITE,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,next->text);
41771f08665SBarry Smith         break;
418b3506946SBarry Smith       case OPTION_STRING:
4191ae3d29cSBarry Smith 	ierr = AMS_Memory_add_field(amem,next->text,next->data,1,AMS_STRING,AMS_WRITE,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,next->text);
4201ae3d29cSBarry Smith         break;
4211ae3d29cSBarry Smith       case OPTION_STRING_ARRAY:
4221ae3d29cSBarry Smith 	ierr = AMS_Memory_add_field(amem,next->text,next->data,next->arraylength,AMS_STRING,AMS_WRITE,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,next->text);
423b3506946SBarry Smith         break;
424b3506946SBarry Smith       case OPTION_LIST:
42571f08665SBarry Smith         {PetscInt  ntext;
42671f08665SBarry Smith         char       ldefault[128];
42771f08665SBarry Smith 	ierr = PetscStrcpy(ldefault,"DEFAULT:");CHKERRQ(ierr);
42871f08665SBarry Smith 	ierr = PetscStrcat(ldefault,next->text);CHKERRQ(ierr);
4291ae3d29cSBarry Smith 	ierr = AMS_Memory_add_field(amem,ldefault,next->data,1,AMS_STRING,AMS_WRITE,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,ldefault);
430065533a5SJed Brown 	ierr = PetscFListGet(next->flist,(const char***)&next->edata,&ntext);CHKERRQ(ierr);
431d5649816SBarry Smith 	ierr = AMS_Memory_add_field(amem,next->text,next->edata,ntext-1,AMS_STRING,AMS_WRITE,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,next->text);
4321ae3d29cSBarry Smith         break;}
4331ae3d29cSBarry Smith       case OPTION_ELIST:
434d5649816SBarry Smith         {PetscInt  ntext = next->nlist;
4351ae3d29cSBarry Smith         char       ldefault[128];
4361ae3d29cSBarry Smith 	ierr = PetscStrcpy(ldefault,"DEFAULT:");CHKERRQ(ierr);
4371ae3d29cSBarry Smith 	ierr = PetscStrcat(ldefault,next->text);CHKERRQ(ierr);
4381ae3d29cSBarry Smith 	ierr = AMS_Memory_add_field(amem,ldefault,next->data,1,AMS_STRING,AMS_WRITE,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,ldefault);
4391ae3d29cSBarry Smith 	ierr = PetscMalloc((ntext+1)*sizeof(char**),&next->edata);CHKERRQ(ierr);
4401ae3d29cSBarry Smith         ierr = PetscMemcpy(next->edata,next->list,ntext*sizeof(char*));CHKERRQ(ierr);
4411ae3d29cSBarry Smith 	ierr = AMS_Memory_add_field(amem,next->text,next->edata,ntext,AMS_STRING,AMS_WRITE,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,next->text);
44271f08665SBarry Smith         break;}
443b3506946SBarry Smith     default:
444b3506946SBarry Smith       break;
445b3506946SBarry Smith     }
446b3506946SBarry Smith     next = next->next;
447b3506946SBarry Smith   }
448b3506946SBarry Smith 
4491ae3d29cSBarry Smith   ierr = AMS_Memory_publish(amem);CHKERRAMS(ierr);
4501ae3d29cSBarry Smith   ierr = AMS_Memory_grant_access(amem);CHKERRAMS(ierr);
451b3506946SBarry Smith   /* wait until accessor has unlocked the memory */
4521ae3d29cSBarry Smith   ierr = AMS_Memory_lock(amem,0);CHKERRAMS(ierr);
4531ae3d29cSBarry Smith   ierr = AMS_Memory_take_access(amem);CHKERRAMS(ierr);
454b3506946SBarry Smith 
455b3506946SBarry Smith   /* reset counter to -2; this updates the screen with the new options for the selected method */
456b3506946SBarry Smith   if (changedmethod) PetscOptionsPublishCount = -2;
457b3506946SBarry Smith 
4581ae3d29cSBarry Smith   ierr = AMS_Memory_grant_access(amem);CHKERRAMS(ierr);
4591ae3d29cSBarry Smith   ierr = AMS_Memory_destroy(amem);CHKERRAMS(ierr);
460b3506946SBarry Smith   PetscFunctionReturn(0);
461b3506946SBarry Smith }
462b3506946SBarry Smith #endif
463b3506946SBarry Smith 
4646356e834SBarry Smith #undef __FUNCT__
46553acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEnd_Private"
46653acd3b1SBarry Smith PetscErrorCode PetscOptionsEnd_Private(void)
46753acd3b1SBarry Smith {
46853acd3b1SBarry Smith   PetscErrorCode ierr;
4696356e834SBarry Smith   PetscOptions   last;
4706356e834SBarry Smith   char           option[256],value[1024],tmp[32];
471330cf3c9SBarry Smith   size_t         j;
47253acd3b1SBarry Smith 
47353acd3b1SBarry Smith   PetscFunctionBegin;
4746356e834SBarry Smith 
4751bc75a8dSBarry Smith   CHKMEMQ;
476aee2cecaSBarry Smith   if (PetscOptionsObject.next) {
477b3506946SBarry Smith     if (!PetscOptionsPublishCount) {
478b3506946SBarry Smith #if defined(PETSC_HAVE_AMS)
479d5649816SBarry Smith       ierr = PetscOptionsAMSInput();CHKERRQ(ierr);
480b3506946SBarry Smith #else
48171f08665SBarry Smith       ierr = PetscOptionsGetFromTextInput();CHKERRQ(ierr);
482b3506946SBarry Smith #endif
483aee2cecaSBarry Smith     }
484aee2cecaSBarry Smith   }
4856356e834SBarry Smith 
486c31cb41cSBarry Smith   ierr = PetscFree(PetscOptionsObject.title);CHKERRQ(ierr);
487c31cb41cSBarry Smith   ierr = PetscFree(PetscOptionsObject.prefix);CHKERRQ(ierr);
4886356e834SBarry Smith 
4896356e834SBarry Smith   /* reset counter to -2; this updates the screen with the new options for the selected method */
4906356e834SBarry Smith   if (PetscOptionsObject.changedmethod) PetscOptionsPublishCount = -2;
49161b37b28SSatish Balay   /* reset alreadyprinted flag */
49261b37b28SSatish Balay   PetscOptionsObject.alreadyprinted = PETSC_FALSE;
4933194b578SJed Brown   if (PetscOptionsObject.object) PetscOptionsObject.object->optionsprinted = PETSC_TRUE;
4943194b578SJed Brown   PetscOptionsObject.object = PETSC_NULL;
4956356e834SBarry Smith 
4966356e834SBarry Smith   while (PetscOptionsObject.next) {
4976356e834SBarry Smith     if (PetscOptionsObject.next->set) {
4986356e834SBarry Smith       if (PetscOptionsObject.prefix) {
4996356e834SBarry Smith         ierr = PetscStrcpy(option,"-");CHKERRQ(ierr);
5006356e834SBarry Smith         ierr = PetscStrcat(option,PetscOptionsObject.prefix);CHKERRQ(ierr);
5016356e834SBarry Smith         ierr = PetscStrcat(option,PetscOptionsObject.next->option+1);CHKERRQ(ierr);
5026356e834SBarry Smith       } else {
5036356e834SBarry Smith         ierr = PetscStrcpy(option,PetscOptionsObject.next->option);CHKERRQ(ierr);
5046356e834SBarry Smith       }
5056356e834SBarry Smith 
5066356e834SBarry Smith       switch (PetscOptionsObject.next->type) {
5076356e834SBarry Smith         case OPTION_HEAD:
5086356e834SBarry Smith           break;
509e26ddf31SBarry Smith         case OPTION_INT_ARRAY:
510e26ddf31SBarry Smith           sprintf(value,"%d",(int)((PetscInt*)PetscOptionsObject.next->data)[0]);
511e26ddf31SBarry Smith           for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
512e26ddf31SBarry Smith             sprintf(tmp,"%d",(int)((PetscInt*)PetscOptionsObject.next->data)[j]);
513e26ddf31SBarry Smith             ierr = PetscStrcat(value,",");CHKERRQ(ierr);
514e26ddf31SBarry Smith             ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
515e26ddf31SBarry Smith           }
516e26ddf31SBarry Smith           break;
5176356e834SBarry Smith         case OPTION_INT:
5187a72a596SBarry Smith           sprintf(value,"%d",(int) *(PetscInt*)PetscOptionsObject.next->data);
5196356e834SBarry Smith           break;
5206356e834SBarry Smith         case OPTION_REAL:
5217a72a596SBarry Smith           sprintf(value,"%g",(double) *(PetscReal*)PetscOptionsObject.next->data);
5226356e834SBarry Smith           break;
5236356e834SBarry Smith         case OPTION_REAL_ARRAY:
5247a72a596SBarry Smith           sprintf(value,"%g",(double)((PetscReal*)PetscOptionsObject.next->data)[0]);
5256356e834SBarry Smith           for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
5267a72a596SBarry Smith             sprintf(tmp,"%g",(double)((PetscReal*)PetscOptionsObject.next->data)[j]);
5276356e834SBarry Smith             ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5286356e834SBarry Smith             ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5296356e834SBarry Smith           }
5306356e834SBarry Smith           break;
5316356e834SBarry Smith         case OPTION_LOGICAL:
53271f08665SBarry Smith           sprintf(value,"%d",*(int*)PetscOptionsObject.next->data);
5336356e834SBarry Smith           break;
5341ae3d29cSBarry Smith       case OPTION_LOGICAL_ARRAY:
535ace3abfcSBarry Smith           sprintf(value,"%d",(int)((PetscBool*)PetscOptionsObject.next->data)[0]);
5361ae3d29cSBarry Smith           for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
537ace3abfcSBarry Smith             sprintf(tmp,"%d",(int)((PetscBool*)PetscOptionsObject.next->data)[j]);
5381ae3d29cSBarry Smith             ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5391ae3d29cSBarry Smith             ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5401ae3d29cSBarry Smith           }
5411ae3d29cSBarry Smith           break;
5426356e834SBarry Smith         case OPTION_LIST:
5431ae3d29cSBarry Smith         case OPTION_ELIST:
54471f08665SBarry Smith           ierr = PetscStrcpy(value,*(char**)PetscOptionsObject.next->data);CHKERRQ(ierr);
5456356e834SBarry Smith           break;
5461ae3d29cSBarry Smith         case OPTION_STRING:
54771f08665SBarry Smith           ierr = PetscStrcpy(value,*(char**)PetscOptionsObject.next->data);CHKERRQ(ierr);
5481ae3d29cSBarry Smith         case OPTION_STRING_ARRAY:
5491ae3d29cSBarry Smith           sprintf(value,"%s",((char**)PetscOptionsObject.next->data)[0]);
5501ae3d29cSBarry Smith           for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
5511ae3d29cSBarry Smith             sprintf(tmp,"%s",((char**)PetscOptionsObject.next->data)[j]);
5521ae3d29cSBarry Smith             ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5531ae3d29cSBarry Smith             ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5541ae3d29cSBarry Smith           }
5556356e834SBarry Smith           break;
5566356e834SBarry Smith       }
5576356e834SBarry Smith       ierr = PetscOptionsSetValue(option,value);CHKERRQ(ierr);
5586356e834SBarry Smith     }
559503cfb0cSBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->text);CHKERRQ(ierr);
560503cfb0cSBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->option);CHKERRQ(ierr);
5616356e834SBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->man);CHKERRQ(ierr);
56205b42c5fSBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->data);CHKERRQ(ierr);
56371f08665SBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->edata);CHKERRQ(ierr);
5646356e834SBarry Smith     last                    = PetscOptionsObject.next;
5656356e834SBarry Smith     PetscOptionsObject.next = PetscOptionsObject.next->next;
5666356e834SBarry Smith     ierr                    = PetscFree(last);CHKERRQ(ierr);
5671bc75a8dSBarry Smith     CHKMEMQ;
5686356e834SBarry Smith   }
5691bc75a8dSBarry Smith   CHKMEMQ;
5706356e834SBarry Smith   PetscOptionsObject.next = 0;
57153acd3b1SBarry Smith   PetscFunctionReturn(0);
57253acd3b1SBarry Smith }
57353acd3b1SBarry Smith 
57453acd3b1SBarry Smith #undef __FUNCT__
57553acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEnum"
57653acd3b1SBarry Smith /*@C
57753acd3b1SBarry Smith    PetscOptionsEnum - Gets the enum value for a particular option in the database.
57853acd3b1SBarry Smith 
5793f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
58053acd3b1SBarry Smith 
58153acd3b1SBarry Smith    Input Parameters:
58253acd3b1SBarry Smith +  opt - option name
58353acd3b1SBarry Smith .  text - short string that describes the option
58453acd3b1SBarry Smith .  man - manual page with additional information on option
58553acd3b1SBarry Smith .  list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null
58653acd3b1SBarry Smith -  defaultv - the default (current) value
58753acd3b1SBarry Smith 
58853acd3b1SBarry Smith    Output Parameter:
58953acd3b1SBarry Smith +  value - the  value to return
59053acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
59153acd3b1SBarry Smith 
59253acd3b1SBarry Smith    Level: beginner
59353acd3b1SBarry Smith 
59453acd3b1SBarry Smith    Concepts: options database
59553acd3b1SBarry Smith 
59653acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
59753acd3b1SBarry Smith 
59853acd3b1SBarry Smith           list is usually something like PCASMTypes or some other predefined list of enum names
59953acd3b1SBarry Smith 
60053acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
601acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
602acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
60353acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
60453acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
605acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
60653acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
60753acd3b1SBarry Smith @*/
6087087cfbeSBarry Smith PetscErrorCode  PetscOptionsEnum(const char opt[],const char text[],const char man[],const char *const*list,PetscEnum defaultv,PetscEnum *value,PetscBool  *set)
60953acd3b1SBarry Smith {
61053acd3b1SBarry Smith   PetscErrorCode ierr;
61153acd3b1SBarry Smith   PetscInt       ntext = 0;
612aa5bb8c0SSatish Balay   PetscInt       tval;
613ace3abfcSBarry Smith   PetscBool      tflg;
61453acd3b1SBarry Smith 
61553acd3b1SBarry Smith   PetscFunctionBegin;
61653acd3b1SBarry Smith   while (list[ntext++]) {
617e32f2f54SBarry Smith     if (ntext > 50) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument appears to be wrong or have more than 50 entries");
61853acd3b1SBarry Smith   }
619e32f2f54SBarry Smith   if (ntext < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument must have at least two entries: typename and type prefix");
62053acd3b1SBarry Smith   ntext -= 3;
621aa5bb8c0SSatish Balay   ierr = PetscOptionsEList(opt,text,man,list,ntext,list[defaultv],&tval,&tflg);CHKERRQ(ierr);
622aa5bb8c0SSatish Balay   /* with PETSC_USE_64BIT_INDICES sizeof(PetscInt) != sizeof(PetscEnum) */
623aa5bb8c0SSatish Balay   if (tflg) *value = (PetscEnum)tval;
624aa5bb8c0SSatish Balay   if (set)  *set   = tflg;
62553acd3b1SBarry Smith   PetscFunctionReturn(0);
62653acd3b1SBarry Smith }
62753acd3b1SBarry Smith 
62853acd3b1SBarry Smith /* -------------------------------------------------------------------------------------------------------------*/
62953acd3b1SBarry Smith #undef __FUNCT__
63053acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsInt"
63153acd3b1SBarry Smith /*@C
63253acd3b1SBarry Smith    PetscOptionsInt - Gets the integer value for a particular option in the database.
63353acd3b1SBarry Smith 
6343f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
63553acd3b1SBarry Smith 
63653acd3b1SBarry Smith    Input Parameters:
63753acd3b1SBarry Smith +  opt - option name
63853acd3b1SBarry Smith .  text - short string that describes the option
63953acd3b1SBarry Smith .  man - manual page with additional information on option
64053acd3b1SBarry Smith -  defaultv - the default (current) value
64153acd3b1SBarry Smith 
64253acd3b1SBarry Smith    Output Parameter:
64353acd3b1SBarry Smith +  value - the integer value to return
64453acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
64553acd3b1SBarry Smith 
64653acd3b1SBarry Smith    Level: beginner
64753acd3b1SBarry Smith 
64853acd3b1SBarry Smith    Concepts: options database^has int
64953acd3b1SBarry Smith 
65053acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
65153acd3b1SBarry Smith 
65253acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
653acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
654acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
65553acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
65653acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
657acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
65853acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
65953acd3b1SBarry Smith @*/
6607087cfbeSBarry Smith PetscErrorCode  PetscOptionsInt(const char opt[],const char text[],const char man[],PetscInt defaultv,PetscInt *value,PetscBool  *set)
66153acd3b1SBarry Smith {
66253acd3b1SBarry Smith   PetscErrorCode ierr;
6636356e834SBarry Smith   PetscOptions   amsopt;
66453acd3b1SBarry Smith 
66553acd3b1SBarry Smith   PetscFunctionBegin;
666b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
6676356e834SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_INT,&amsopt);CHKERRQ(ierr);
6686356e834SBarry Smith     ierr = PetscMalloc(sizeof(PetscInt),&amsopt->data);CHKERRQ(ierr);
6696356e834SBarry Smith     *(PetscInt*)amsopt->data = defaultv;
670af6d86caSBarry Smith   }
67153acd3b1SBarry Smith   ierr = PetscOptionsGetInt(PetscOptionsObject.prefix,opt,value,set);CHKERRQ(ierr);
67261b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
6732aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%d>: %s (%s)\n",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,defaultv,text,ManSection(man));CHKERRQ(ierr);
67453acd3b1SBarry Smith   }
67553acd3b1SBarry Smith   PetscFunctionReturn(0);
67653acd3b1SBarry Smith }
67753acd3b1SBarry Smith 
67853acd3b1SBarry Smith #undef __FUNCT__
67953acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsString"
68053acd3b1SBarry Smith /*@C
68153acd3b1SBarry Smith    PetscOptionsString - Gets the string value for a particular option in the database.
68253acd3b1SBarry Smith 
6833f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
68453acd3b1SBarry Smith 
68553acd3b1SBarry Smith    Input Parameters:
68653acd3b1SBarry Smith +  opt - option name
68753acd3b1SBarry Smith .  text - short string that describes the option
68853acd3b1SBarry Smith .  man - manual page with additional information on option
689bcbf2dc5SJed Brown .  defaultv - the default (current) value
690bcbf2dc5SJed Brown -  len - length of the result string including null terminator
69153acd3b1SBarry Smith 
69253acd3b1SBarry Smith    Output Parameter:
69353acd3b1SBarry Smith +  value - the value to return
69453acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
69553acd3b1SBarry Smith 
69653acd3b1SBarry Smith    Level: beginner
69753acd3b1SBarry Smith 
69853acd3b1SBarry Smith    Concepts: options database^has int
69953acd3b1SBarry Smith 
70053acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
70153acd3b1SBarry Smith 
7027fccdfe4SBarry 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).
7037fccdfe4SBarry Smith 
70453acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
705acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
706acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsReal(), PetscOptionsBool(),
70753acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
70853acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
709acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
71053acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
71153acd3b1SBarry Smith @*/
7127087cfbeSBarry Smith PetscErrorCode  PetscOptionsString(const char opt[],const char text[],const char man[],const char defaultv[],char value[],size_t len,PetscBool  *set)
71353acd3b1SBarry Smith {
71453acd3b1SBarry Smith   PetscErrorCode ierr;
715aee2cecaSBarry Smith   PetscOptions   amsopt;
71653acd3b1SBarry Smith 
71753acd3b1SBarry Smith   PetscFunctionBegin;
718b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
719aee2cecaSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
72071f08665SBarry Smith     ierr = PetscMalloc(sizeof(char*),&amsopt->data);CHKERRQ(ierr);
72171f08665SBarry Smith     *(const char**)amsopt->data = defaultv;
722af6d86caSBarry Smith   }
72353acd3b1SBarry Smith   ierr = PetscOptionsGetString(PetscOptionsObject.prefix,opt,value,len,set);CHKERRQ(ierr);
72461b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
7252aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%s>: %s (%s)\n",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,defaultv,text,ManSection(man));CHKERRQ(ierr);
72653acd3b1SBarry Smith   }
72753acd3b1SBarry Smith   PetscFunctionReturn(0);
72853acd3b1SBarry Smith }
72953acd3b1SBarry Smith 
73053acd3b1SBarry Smith #undef __FUNCT__
73153acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsReal"
73253acd3b1SBarry Smith /*@C
73353acd3b1SBarry Smith    PetscOptionsReal - Gets the PetscReal value for a particular option in the database.
73453acd3b1SBarry Smith 
7353f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
73653acd3b1SBarry Smith 
73753acd3b1SBarry Smith    Input Parameters:
73853acd3b1SBarry Smith +  opt - option name
73953acd3b1SBarry Smith .  text - short string that describes the option
74053acd3b1SBarry Smith .  man - manual page with additional information on option
74153acd3b1SBarry Smith -  defaultv - the default (current) value
74253acd3b1SBarry Smith 
74353acd3b1SBarry Smith    Output Parameter:
74453acd3b1SBarry Smith +  value - the value to return
74553acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
74653acd3b1SBarry Smith 
74753acd3b1SBarry Smith    Level: beginner
74853acd3b1SBarry Smith 
74953acd3b1SBarry Smith    Concepts: options database^has int
75053acd3b1SBarry Smith 
75153acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
75253acd3b1SBarry Smith 
75353acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
754acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
755acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
75653acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
75753acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
758acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
75953acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
76053acd3b1SBarry Smith @*/
7617087cfbeSBarry Smith PetscErrorCode  PetscOptionsReal(const char opt[],const char text[],const char man[],PetscReal defaultv,PetscReal *value,PetscBool  *set)
76253acd3b1SBarry Smith {
76353acd3b1SBarry Smith   PetscErrorCode ierr;
764538aa990SBarry Smith   PetscOptions   amsopt;
76553acd3b1SBarry Smith 
76653acd3b1SBarry Smith   PetscFunctionBegin;
767b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
768538aa990SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_REAL,&amsopt);CHKERRQ(ierr);
769538aa990SBarry Smith     ierr = PetscMalloc(sizeof(PetscReal),&amsopt->data);CHKERRQ(ierr);
770538aa990SBarry Smith     *(PetscReal*)amsopt->data = defaultv;
771538aa990SBarry Smith   }
77253acd3b1SBarry Smith   ierr = PetscOptionsGetReal(PetscOptionsObject.prefix,opt,value,set);CHKERRQ(ierr);
77361b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
7742aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%G>: %s (%s)\n",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,defaultv,text,ManSection(man));CHKERRQ(ierr);
77553acd3b1SBarry Smith   }
77653acd3b1SBarry Smith   PetscFunctionReturn(0);
77753acd3b1SBarry Smith }
77853acd3b1SBarry Smith 
77953acd3b1SBarry Smith #undef __FUNCT__
78053acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsScalar"
78153acd3b1SBarry Smith /*@C
78253acd3b1SBarry Smith    PetscOptionsScalar - Gets the scalar value for a particular option in the database.
78353acd3b1SBarry Smith 
7843f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
78553acd3b1SBarry Smith 
78653acd3b1SBarry Smith    Input Parameters:
78753acd3b1SBarry Smith +  opt - option name
78853acd3b1SBarry Smith .  text - short string that describes the option
78953acd3b1SBarry Smith .  man - manual page with additional information on option
79053acd3b1SBarry Smith -  defaultv - the default (current) value
79153acd3b1SBarry Smith 
79253acd3b1SBarry Smith    Output Parameter:
79353acd3b1SBarry Smith +  value - the value to return
79453acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
79553acd3b1SBarry Smith 
79653acd3b1SBarry Smith    Level: beginner
79753acd3b1SBarry Smith 
79853acd3b1SBarry Smith    Concepts: options database^has int
79953acd3b1SBarry Smith 
80053acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
80153acd3b1SBarry Smith 
80253acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
803acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
804acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
80553acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
80653acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
807acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
80853acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
80953acd3b1SBarry Smith @*/
8107087cfbeSBarry Smith PetscErrorCode  PetscOptionsScalar(const char opt[],const char text[],const char man[],PetscScalar defaultv,PetscScalar *value,PetscBool  *set)
81153acd3b1SBarry Smith {
81253acd3b1SBarry Smith   PetscErrorCode ierr;
81353acd3b1SBarry Smith 
81453acd3b1SBarry Smith   PetscFunctionBegin;
81553acd3b1SBarry Smith #if !defined(PETSC_USE_COMPLEX)
81653acd3b1SBarry Smith   ierr = PetscOptionsReal(opt,text,man,defaultv,value,set);CHKERRQ(ierr);
81753acd3b1SBarry Smith #else
81853acd3b1SBarry Smith   ierr = PetscOptionsGetScalar(PetscOptionsObject.prefix,opt,value,set);CHKERRQ(ierr);
81953acd3b1SBarry Smith #endif
82053acd3b1SBarry Smith   PetscFunctionReturn(0);
82153acd3b1SBarry Smith }
82253acd3b1SBarry Smith 
82353acd3b1SBarry Smith #undef __FUNCT__
82453acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsName"
82553acd3b1SBarry Smith /*@C
82690d69ab7SBarry 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
82790d69ab7SBarry Smith                       its value is set to false.
82853acd3b1SBarry Smith 
8293f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
83053acd3b1SBarry Smith 
83153acd3b1SBarry Smith    Input Parameters:
83253acd3b1SBarry Smith +  opt - option name
83353acd3b1SBarry Smith .  text - short string that describes the option
83453acd3b1SBarry Smith -  man - manual page with additional information on option
83553acd3b1SBarry Smith 
83653acd3b1SBarry Smith    Output Parameter:
83753acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
83853acd3b1SBarry Smith 
83953acd3b1SBarry Smith    Level: beginner
84053acd3b1SBarry Smith 
84153acd3b1SBarry Smith    Concepts: options database^has int
84253acd3b1SBarry Smith 
84353acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
84453acd3b1SBarry Smith 
84553acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
846acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
847acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
84853acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
84953acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
850acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
85153acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
85253acd3b1SBarry Smith @*/
8537087cfbeSBarry Smith PetscErrorCode  PetscOptionsName(const char opt[],const char text[],const char man[],PetscBool  *flg)
85453acd3b1SBarry Smith {
85553acd3b1SBarry Smith   PetscErrorCode ierr;
8561ae3d29cSBarry Smith   PetscOptions   amsopt;
85753acd3b1SBarry Smith 
85853acd3b1SBarry Smith   PetscFunctionBegin;
8591ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
8601ae3d29cSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_LOGICAL,&amsopt);CHKERRQ(ierr);
861ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
862ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
8631ae3d29cSBarry Smith   }
86453acd3b1SBarry Smith   ierr = PetscOptionsHasName(PetscOptionsObject.prefix,opt,flg);CHKERRQ(ierr);
86561b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
8662aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s: %s (%s)\n",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,text,ManSection(man));CHKERRQ(ierr);
86753acd3b1SBarry Smith   }
86853acd3b1SBarry Smith   PetscFunctionReturn(0);
86953acd3b1SBarry Smith }
87053acd3b1SBarry Smith 
87153acd3b1SBarry Smith #undef __FUNCT__
87253acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsList"
87353acd3b1SBarry Smith /*@C
87453acd3b1SBarry Smith      PetscOptionsList - Puts a list of option values that a single one may be selected from
87553acd3b1SBarry Smith 
8763f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
87753acd3b1SBarry Smith 
87853acd3b1SBarry Smith    Input Parameters:
87953acd3b1SBarry Smith +  opt - option name
88053acd3b1SBarry Smith .  text - short string that describes the option
88153acd3b1SBarry Smith .  man - manual page with additional information on option
88253acd3b1SBarry Smith .  list - the possible choices
8833cc1e11dSBarry Smith .  defaultv - the default (current) value
8843cc1e11dSBarry Smith -  len - the length of the character array value
88553acd3b1SBarry Smith 
88653acd3b1SBarry Smith    Output Parameter:
88753acd3b1SBarry Smith +  value - the value to return
88853acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
88953acd3b1SBarry Smith 
89053acd3b1SBarry Smith    Level: intermediate
89153acd3b1SBarry Smith 
89253acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
89353acd3b1SBarry Smith 
89453acd3b1SBarry Smith    See PetscOptionsEList() for when the choices are given in a string array
89553acd3b1SBarry Smith 
89653acd3b1SBarry Smith    To get a listing of all currently specified options,
89788c29154SBarry Smith     see PetscOptionsView() or PetscOptionsGetAll()
89853acd3b1SBarry Smith 
89953acd3b1SBarry Smith    Concepts: options database^list
90053acd3b1SBarry Smith 
90153acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
902acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
90353acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
90453acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
905acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
906171400e9SBarry Smith           PetscOptionsList(), PetscOptionsEList(), PetscOptionsEnum()
90753acd3b1SBarry Smith @*/
9087087cfbeSBarry Smith PetscErrorCode  PetscOptionsList(const char opt[],const char ltext[],const char man[],PetscFList list,const char defaultv[],char value[],size_t len,PetscBool  *set)
90953acd3b1SBarry Smith {
91053acd3b1SBarry Smith   PetscErrorCode ierr;
9113cc1e11dSBarry Smith   PetscOptions   amsopt;
91253acd3b1SBarry Smith 
91353acd3b1SBarry Smith   PetscFunctionBegin;
914b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
9153cc1e11dSBarry Smith     ierr = PetscOptionsCreate_Private(opt,ltext,man,OPTION_LIST,&amsopt);CHKERRQ(ierr);
91671f08665SBarry Smith     ierr = PetscMalloc(sizeof(char*),&amsopt->data);CHKERRQ(ierr);
91771f08665SBarry Smith     *(const char**)amsopt->data = defaultv;
9183cc1e11dSBarry Smith     amsopt->flist = list;
9193cc1e11dSBarry Smith   }
92053acd3b1SBarry Smith   ierr = PetscOptionsGetString(PetscOptionsObject.prefix,opt,value,len,set);CHKERRQ(ierr);
92161b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
9223cc1e11dSBarry Smith     ierr = PetscFListPrintTypes(PetscOptionsObject.comm,stdout,PetscOptionsObject.prefix,opt,ltext,man,list,defaultv);CHKERRQ(ierr);CHKERRQ(ierr);
92353acd3b1SBarry Smith   }
92453acd3b1SBarry Smith   PetscFunctionReturn(0);
92553acd3b1SBarry Smith }
92653acd3b1SBarry Smith 
92753acd3b1SBarry Smith #undef __FUNCT__
92853acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEList"
92953acd3b1SBarry Smith /*@C
93053acd3b1SBarry Smith      PetscOptionsEList - Puts a list of option values that a single one may be selected from
93153acd3b1SBarry Smith 
9323f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
93353acd3b1SBarry Smith 
93453acd3b1SBarry Smith    Input Parameters:
93553acd3b1SBarry Smith +  opt - option name
93653acd3b1SBarry Smith .  ltext - short string that describes the option
93753acd3b1SBarry Smith .  man - manual page with additional information on option
93853acd3b1SBarry Smith .  list - the possible choices
93953acd3b1SBarry Smith .  ntext - number of choices
94053acd3b1SBarry Smith -  defaultv - the default (current) value
94153acd3b1SBarry Smith 
94253acd3b1SBarry Smith    Output Parameter:
94353acd3b1SBarry Smith +  value - the index of the value to return
94453acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
94553acd3b1SBarry Smith 
94653acd3b1SBarry Smith    Level: intermediate
94753acd3b1SBarry Smith 
94853acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
94953acd3b1SBarry Smith 
95053acd3b1SBarry Smith    See PetscOptionsList() for when the choices are given in a PetscFList()
95153acd3b1SBarry Smith 
95253acd3b1SBarry Smith    Concepts: options database^list
95353acd3b1SBarry Smith 
95453acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
955acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
95653acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
95753acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
958acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
959171400e9SBarry Smith           PetscOptionsList(), PetscOptionsEnum()
96053acd3b1SBarry Smith @*/
9617087cfbeSBarry 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)
96253acd3b1SBarry Smith {
96353acd3b1SBarry Smith   PetscErrorCode ierr;
96453acd3b1SBarry Smith   PetscInt       i;
9651ae3d29cSBarry Smith   PetscOptions   amsopt;
96653acd3b1SBarry Smith 
96753acd3b1SBarry Smith   PetscFunctionBegin;
9681ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
969d5649816SBarry Smith     ierr = PetscOptionsCreate_Private(opt,ltext,man,OPTION_ELIST,&amsopt);CHKERRQ(ierr);
9701ae3d29cSBarry Smith     ierr = PetscMalloc(sizeof(char*),&amsopt->data);CHKERRQ(ierr);
9711ae3d29cSBarry Smith     *(const char**)amsopt->data = defaultv;
9721ae3d29cSBarry Smith     amsopt->list  = list;
9731ae3d29cSBarry Smith     amsopt->nlist = ntext;
9741ae3d29cSBarry Smith   }
97553acd3b1SBarry Smith   ierr = PetscOptionsGetEList(PetscOptionsObject.prefix,opt,list,ntext,value,set);CHKERRQ(ierr);
97661b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
97753acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%s> (choose one of)",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,defaultv);CHKERRQ(ierr);
97853acd3b1SBarry Smith     for (i=0; i<ntext; i++){
97953acd3b1SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm," %s",list[i]);CHKERRQ(ierr);
98053acd3b1SBarry Smith     }
9812aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm," (%s)\n",ManSection(man));CHKERRQ(ierr);
98253acd3b1SBarry Smith   }
98353acd3b1SBarry Smith   PetscFunctionReturn(0);
98453acd3b1SBarry Smith }
98553acd3b1SBarry Smith 
98653acd3b1SBarry Smith #undef __FUNCT__
987acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolGroupBegin"
98853acd3b1SBarry Smith /*@C
989acfcf0e5SJed Brown      PetscOptionsBoolGroupBegin - First in a series of logical queries on the options database for
990d5649816SBarry Smith        which at most a single value can be true.
99153acd3b1SBarry Smith 
9923f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
99353acd3b1SBarry Smith 
99453acd3b1SBarry Smith    Input Parameters:
99553acd3b1SBarry Smith +  opt - option name
99653acd3b1SBarry Smith .  text - short string that describes the option
99753acd3b1SBarry Smith -  man - manual page with additional information on option
99853acd3b1SBarry Smith 
99953acd3b1SBarry Smith    Output Parameter:
100053acd3b1SBarry Smith .  flg - whether that option was set or not
100153acd3b1SBarry Smith 
100253acd3b1SBarry Smith    Level: intermediate
100353acd3b1SBarry Smith 
100453acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
100553acd3b1SBarry Smith 
1006acfcf0e5SJed Brown    Must be followed by 0 or more PetscOptionsBoolGroup()s and PetscOptionsBoolGroupEnd()
100753acd3b1SBarry Smith 
100853acd3b1SBarry Smith     Concepts: options database^logical group
100953acd3b1SBarry Smith 
101053acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1011acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
101253acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
101353acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1014acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
101553acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
101653acd3b1SBarry Smith @*/
10177087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolGroupBegin(const char opt[],const char text[],const char man[],PetscBool  *flg)
101853acd3b1SBarry Smith {
101953acd3b1SBarry Smith   PetscErrorCode ierr;
10201ae3d29cSBarry Smith   PetscOptions   amsopt;
102153acd3b1SBarry Smith 
102253acd3b1SBarry Smith   PetscFunctionBegin;
10231ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
10241ae3d29cSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_LOGICAL,&amsopt);CHKERRQ(ierr);
1025ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1026ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
10271ae3d29cSBarry Smith   }
102868b16fdaSBarry Smith   *flg = PETSC_FALSE;
1029acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,PETSC_NULL);CHKERRQ(ierr);
103061b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
103153acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  Pick at most one of -------------\n");CHKERRQ(ierr);
10322aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"    -%s%s: %s (%s)\n",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,text,ManSection(man));CHKERRQ(ierr);
103353acd3b1SBarry Smith   }
103453acd3b1SBarry Smith   PetscFunctionReturn(0);
103553acd3b1SBarry Smith }
103653acd3b1SBarry Smith 
103753acd3b1SBarry Smith #undef __FUNCT__
1038acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolGroup"
103953acd3b1SBarry Smith /*@C
1040acfcf0e5SJed Brown      PetscOptionsBoolGroup - One in a series of logical queries on the options database for
1041d5649816SBarry Smith        which at most a single value can be true.
104253acd3b1SBarry Smith 
10433f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
104453acd3b1SBarry Smith 
104553acd3b1SBarry Smith    Input Parameters:
104653acd3b1SBarry Smith +  opt - option name
104753acd3b1SBarry Smith .  text - short string that describes the option
104853acd3b1SBarry Smith -  man - manual page with additional information on option
104953acd3b1SBarry Smith 
105053acd3b1SBarry Smith    Output Parameter:
105153acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
105253acd3b1SBarry Smith 
105353acd3b1SBarry Smith    Level: intermediate
105453acd3b1SBarry Smith 
105553acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
105653acd3b1SBarry Smith 
1057acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin() and preceded a PetscOptionsBoolGroupEnd()
105853acd3b1SBarry Smith 
105953acd3b1SBarry Smith     Concepts: options database^logical group
106053acd3b1SBarry Smith 
106153acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1062acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
106353acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
106453acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1065acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
106653acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
106753acd3b1SBarry Smith @*/
10687087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolGroup(const char opt[],const char text[],const char man[],PetscBool  *flg)
106953acd3b1SBarry Smith {
107053acd3b1SBarry Smith   PetscErrorCode ierr;
10711ae3d29cSBarry Smith   PetscOptions   amsopt;
107253acd3b1SBarry Smith 
107353acd3b1SBarry Smith   PetscFunctionBegin;
10741ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
10751ae3d29cSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_LOGICAL,&amsopt);CHKERRQ(ierr);
1076ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1077ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
10781ae3d29cSBarry Smith   }
107917326d04SJed Brown   *flg = PETSC_FALSE;
1080acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,PETSC_NULL);CHKERRQ(ierr);
108161b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
10822aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"    -%s%s: %s (%s)\n",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,text,ManSection(man));CHKERRQ(ierr);
108353acd3b1SBarry Smith   }
108453acd3b1SBarry Smith   PetscFunctionReturn(0);
108553acd3b1SBarry Smith }
108653acd3b1SBarry Smith 
108753acd3b1SBarry Smith #undef __FUNCT__
1088acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolGroupEnd"
108953acd3b1SBarry Smith /*@C
1090acfcf0e5SJed Brown      PetscOptionsBoolGroupEnd - Last in a series of logical queries on the options database for
1091d5649816SBarry Smith        which at most a single value can be true.
109253acd3b1SBarry Smith 
10933f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
109453acd3b1SBarry Smith 
109553acd3b1SBarry Smith    Input Parameters:
109653acd3b1SBarry Smith +  opt - option name
109753acd3b1SBarry Smith .  text - short string that describes the option
109853acd3b1SBarry Smith -  man - manual page with additional information on option
109953acd3b1SBarry Smith 
110053acd3b1SBarry Smith    Output Parameter:
110153acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
110253acd3b1SBarry Smith 
110353acd3b1SBarry Smith    Level: intermediate
110453acd3b1SBarry Smith 
110553acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
110653acd3b1SBarry Smith 
1107acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin()
110853acd3b1SBarry Smith 
110953acd3b1SBarry Smith     Concepts: options database^logical group
111053acd3b1SBarry Smith 
111153acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1112acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
111353acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
111453acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1115acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
111653acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
111753acd3b1SBarry Smith @*/
11187087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolGroupEnd(const char opt[],const char text[],const char man[],PetscBool  *flg)
111953acd3b1SBarry Smith {
112053acd3b1SBarry Smith   PetscErrorCode ierr;
11211ae3d29cSBarry Smith   PetscOptions   amsopt;
112253acd3b1SBarry Smith 
112353acd3b1SBarry Smith   PetscFunctionBegin;
11241ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
11251ae3d29cSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_LOGICAL,&amsopt);CHKERRQ(ierr);
1126ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1127ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
11281ae3d29cSBarry Smith   }
112917326d04SJed Brown   *flg = PETSC_FALSE;
1130acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,PETSC_NULL);CHKERRQ(ierr);
113161b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
11322aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"    -%s%s: %s (%s)\n",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,text,ManSection(man));CHKERRQ(ierr);
113353acd3b1SBarry Smith   }
113453acd3b1SBarry Smith   PetscFunctionReturn(0);
113553acd3b1SBarry Smith }
113653acd3b1SBarry Smith 
113753acd3b1SBarry Smith #undef __FUNCT__
1138acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBool"
113953acd3b1SBarry Smith /*@C
1140acfcf0e5SJed Brown    PetscOptionsBool - Determines if a particular option is in the database with a true or false
114153acd3b1SBarry Smith 
11423f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
114353acd3b1SBarry Smith 
114453acd3b1SBarry Smith    Input Parameters:
114553acd3b1SBarry Smith +  opt - option name
114653acd3b1SBarry Smith .  text - short string that describes the option
114753acd3b1SBarry Smith -  man - manual page with additional information on option
114853acd3b1SBarry Smith 
114953acd3b1SBarry Smith    Output Parameter:
115053acd3b1SBarry Smith .  flg - PETSC_TRUE or PETSC_FALSE
115153acd3b1SBarry Smith .  set - PETSC_TRUE if found, else PETSC_FALSE
115253acd3b1SBarry Smith 
115353acd3b1SBarry Smith    Level: beginner
115453acd3b1SBarry Smith 
115553acd3b1SBarry Smith    Concepts: options database^logical
115653acd3b1SBarry Smith 
115753acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
115853acd3b1SBarry Smith 
115953acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
1160acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
1161acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
116253acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
116353acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1164acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
116553acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
116653acd3b1SBarry Smith @*/
11677087cfbeSBarry Smith PetscErrorCode  PetscOptionsBool(const char opt[],const char text[],const char man[],PetscBool  deflt,PetscBool  *flg,PetscBool  *set)
116853acd3b1SBarry Smith {
116953acd3b1SBarry Smith   PetscErrorCode ierr;
1170ace3abfcSBarry Smith   PetscBool      iset;
1171aee2cecaSBarry Smith   PetscOptions   amsopt;
117253acd3b1SBarry Smith 
117353acd3b1SBarry Smith   PetscFunctionBegin;
1174b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
1175aee2cecaSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_LOGICAL,&amsopt);CHKERRQ(ierr);
1176ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1177ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = deflt;
1178af6d86caSBarry Smith   }
1179acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,&iset);CHKERRQ(ierr);
118053acd3b1SBarry Smith   if (!iset) {
118153acd3b1SBarry Smith     if (flg) *flg = deflt;
118253acd3b1SBarry Smith   }
118353acd3b1SBarry Smith   if (set) *set = iset;
118461b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
1185ace3abfcSBarry Smith     const char *v = PetscBools[deflt];
11862aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s: <%s> %s (%s)\n",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,v,text,ManSection(man));CHKERRQ(ierr);
118753acd3b1SBarry Smith   }
118853acd3b1SBarry Smith   PetscFunctionReturn(0);
118953acd3b1SBarry Smith }
119053acd3b1SBarry Smith 
119153acd3b1SBarry Smith #undef __FUNCT__
119253acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsRealArray"
119353acd3b1SBarry Smith /*@C
119453acd3b1SBarry Smith    PetscOptionsRealArray - Gets an array of double values for a particular
119553acd3b1SBarry Smith    option in the database. The values must be separated with commas with
119653acd3b1SBarry Smith    no intervening spaces.
119753acd3b1SBarry Smith 
11983f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
119953acd3b1SBarry Smith 
120053acd3b1SBarry Smith    Input Parameters:
120153acd3b1SBarry Smith +  opt - the option one is seeking
120253acd3b1SBarry Smith .  text - short string describing option
120353acd3b1SBarry Smith .  man - manual page for option
120453acd3b1SBarry Smith -  nmax - maximum number of values
120553acd3b1SBarry Smith 
120653acd3b1SBarry Smith    Output Parameter:
120753acd3b1SBarry Smith +  value - location to copy values
120853acd3b1SBarry Smith .  nmax - actual number of values found
120953acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
121053acd3b1SBarry Smith 
121153acd3b1SBarry Smith    Level: beginner
121253acd3b1SBarry Smith 
121353acd3b1SBarry Smith    Notes:
121453acd3b1SBarry Smith    The user should pass in an array of doubles
121553acd3b1SBarry Smith 
121653acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
121753acd3b1SBarry Smith 
121853acd3b1SBarry Smith    Concepts: options database^array of strings
121953acd3b1SBarry Smith 
122053acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1221acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
122253acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
122353acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1224acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
122553acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
122653acd3b1SBarry Smith @*/
12277087cfbeSBarry Smith PetscErrorCode  PetscOptionsRealArray(const char opt[],const char text[],const char man[],PetscReal value[],PetscInt *n,PetscBool  *set)
122853acd3b1SBarry Smith {
122953acd3b1SBarry Smith   PetscErrorCode ierr;
123053acd3b1SBarry Smith   PetscInt       i;
1231e26ddf31SBarry Smith   PetscOptions   amsopt;
123253acd3b1SBarry Smith 
123353acd3b1SBarry Smith   PetscFunctionBegin;
1234b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
1235e26ddf31SBarry Smith     PetscReal *vals;
1236e26ddf31SBarry Smith 
1237e26ddf31SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_REAL_ARRAY,&amsopt);CHKERRQ(ierr);
1238e26ddf31SBarry Smith     ierr = PetscMalloc((*n)*sizeof(PetscReal),&amsopt->data);CHKERRQ(ierr);
1239e26ddf31SBarry Smith     vals = (PetscReal*)amsopt->data;
1240e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1241e26ddf31SBarry Smith     amsopt->arraylength = *n;
1242e26ddf31SBarry Smith   }
124353acd3b1SBarry Smith   ierr = PetscOptionsGetRealArray(PetscOptionsObject.prefix,opt,value,n,set);CHKERRQ(ierr);
124461b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
1245a83599f4SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%G",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,value[0]);CHKERRQ(ierr);
124653acd3b1SBarry Smith     for (i=1; i<*n; i++) {
1247a83599f4SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,",%G",value[i]);CHKERRQ(ierr);
124853acd3b1SBarry Smith     }
12492aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
125053acd3b1SBarry Smith   }
125153acd3b1SBarry Smith   PetscFunctionReturn(0);
125253acd3b1SBarry Smith }
125353acd3b1SBarry Smith 
125453acd3b1SBarry Smith 
125553acd3b1SBarry Smith #undef __FUNCT__
125653acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsIntArray"
125753acd3b1SBarry Smith /*@C
125853acd3b1SBarry Smith    PetscOptionsIntArray - Gets an array of integers for a particular
1259b32a342fSShri Abhyankar    option in the database.
126053acd3b1SBarry Smith 
12613f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
126253acd3b1SBarry Smith 
126353acd3b1SBarry Smith    Input Parameters:
126453acd3b1SBarry Smith +  opt - the option one is seeking
126553acd3b1SBarry Smith .  text - short string describing option
126653acd3b1SBarry Smith .  man - manual page for option
1267f8a50e2bSBarry Smith -  n - maximum number of values
126853acd3b1SBarry Smith 
126953acd3b1SBarry Smith    Output Parameter:
127053acd3b1SBarry Smith +  value - location to copy values
1271f8a50e2bSBarry Smith .  n - actual number of values found
127253acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
127353acd3b1SBarry Smith 
127453acd3b1SBarry Smith    Level: beginner
127553acd3b1SBarry Smith 
127653acd3b1SBarry Smith    Notes:
1277b32a342fSShri Abhyankar    The array can be passed as
1278b32a342fSShri Abhyankar    a comma seperated list:                                 0,1,2,3,4,5,6,7
12790fd488f5SShri Abhyankar    a range (start-end+1):                                  0-8
12800fd488f5SShri Abhyankar    a range with given increment (start-end+1:inc):         0-7:2
12810fd488f5SShri Abhyankar    a combination of values and ranges seperated by commas: 0,1-8,8-15:2
1282b32a342fSShri Abhyankar 
1283b32a342fSShri Abhyankar    There must be no intervening spaces between the values.
128453acd3b1SBarry Smith 
128553acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
128653acd3b1SBarry Smith 
1287b32a342fSShri Abhyankar    Concepts: options database^array of ints
128853acd3b1SBarry Smith 
128953acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1290acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
129153acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
129253acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1293acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
129453acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList(), PetscOptionsRealArray()
129553acd3b1SBarry Smith @*/
12967087cfbeSBarry Smith PetscErrorCode  PetscOptionsIntArray(const char opt[],const char text[],const char man[],PetscInt value[],PetscInt *n,PetscBool  *set)
129753acd3b1SBarry Smith {
129853acd3b1SBarry Smith   PetscErrorCode ierr;
129953acd3b1SBarry Smith   PetscInt       i;
1300e26ddf31SBarry Smith   PetscOptions   amsopt;
130153acd3b1SBarry Smith 
130253acd3b1SBarry Smith   PetscFunctionBegin;
1303b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
1304e26ddf31SBarry Smith     PetscInt *vals;
1305e26ddf31SBarry Smith 
1306e26ddf31SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_INT_ARRAY,&amsopt);CHKERRQ(ierr);
1307e26ddf31SBarry Smith     ierr = PetscMalloc((*n)*sizeof(PetscInt),&amsopt->data);CHKERRQ(ierr);
1308e26ddf31SBarry Smith     vals = (PetscInt*)amsopt->data;
1309e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1310e26ddf31SBarry Smith     amsopt->arraylength = *n;
1311e26ddf31SBarry Smith   }
131253acd3b1SBarry Smith   ierr = PetscOptionsGetIntArray(PetscOptionsObject.prefix,opt,value,n,set);CHKERRQ(ierr);
131361b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
131453acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%d",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,value[0]);CHKERRQ(ierr);
131553acd3b1SBarry Smith     for (i=1; i<*n; i++) {
131653acd3b1SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,",%d",value[i]);CHKERRQ(ierr);
131753acd3b1SBarry Smith     }
13182aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
131953acd3b1SBarry Smith   }
132053acd3b1SBarry Smith   PetscFunctionReturn(0);
132153acd3b1SBarry Smith }
132253acd3b1SBarry Smith 
132353acd3b1SBarry Smith #undef __FUNCT__
132453acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsStringArray"
132553acd3b1SBarry Smith /*@C
132653acd3b1SBarry Smith    PetscOptionsStringArray - Gets an array of string values for a particular
132753acd3b1SBarry Smith    option in the database. The values must be separated with commas with
132853acd3b1SBarry Smith    no intervening spaces.
132953acd3b1SBarry Smith 
13303f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
133153acd3b1SBarry Smith 
133253acd3b1SBarry Smith    Input Parameters:
133353acd3b1SBarry Smith +  opt - the option one is seeking
133453acd3b1SBarry Smith .  text - short string describing option
133553acd3b1SBarry Smith .  man - manual page for option
133653acd3b1SBarry Smith -  nmax - maximum number of strings
133753acd3b1SBarry Smith 
133853acd3b1SBarry Smith    Output Parameter:
133953acd3b1SBarry Smith +  value - location to copy strings
134053acd3b1SBarry Smith .  nmax - actual number of strings found
134153acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
134253acd3b1SBarry Smith 
134353acd3b1SBarry Smith    Level: beginner
134453acd3b1SBarry Smith 
134553acd3b1SBarry Smith    Notes:
134653acd3b1SBarry Smith    The user should pass in an array of pointers to char, to hold all the
134753acd3b1SBarry Smith    strings returned by this function.
134853acd3b1SBarry Smith 
134953acd3b1SBarry Smith    The user is responsible for deallocating the strings that are
135053acd3b1SBarry Smith    returned. The Fortran interface for this routine is not supported.
135153acd3b1SBarry Smith 
135253acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
135353acd3b1SBarry Smith 
135453acd3b1SBarry Smith    Concepts: options database^array of strings
135553acd3b1SBarry Smith 
135653acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1357acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
135853acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
135953acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1360acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
136153acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
136253acd3b1SBarry Smith @*/
13637087cfbeSBarry Smith PetscErrorCode  PetscOptionsStringArray(const char opt[],const char text[],const char man[],char *value[],PetscInt *nmax,PetscBool  *set)
136453acd3b1SBarry Smith {
136553acd3b1SBarry Smith   PetscErrorCode ierr;
13661ae3d29cSBarry Smith   PetscOptions   amsopt;
136753acd3b1SBarry Smith 
136853acd3b1SBarry Smith   PetscFunctionBegin;
13691ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
13701ae3d29cSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_STRING_ARRAY,&amsopt);CHKERRQ(ierr);
13711ae3d29cSBarry Smith     ierr = PetscMalloc((*nmax)*sizeof(char*),&amsopt->data);CHKERRQ(ierr);
13721ae3d29cSBarry Smith     amsopt->arraylength = *nmax;
13731ae3d29cSBarry Smith   }
137453acd3b1SBarry Smith   ierr = PetscOptionsGetStringArray(PetscOptionsObject.prefix,opt,value,nmax,set);CHKERRQ(ierr);
137561b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
13762aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <string1,string2,...>: %s (%s)\n",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,text,ManSection(man));CHKERRQ(ierr);
137753acd3b1SBarry Smith   }
137853acd3b1SBarry Smith   PetscFunctionReturn(0);
137953acd3b1SBarry Smith }
138053acd3b1SBarry Smith 
1381e2446a98SMatthew Knepley #undef __FUNCT__
1382acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolArray"
1383e2446a98SMatthew Knepley /*@C
1384acfcf0e5SJed Brown    PetscOptionsBoolArray - Gets an array of logical values (true or false) for a particular
1385e2446a98SMatthew Knepley    option in the database. The values must be separated with commas with
1386e2446a98SMatthew Knepley    no intervening spaces.
1387e2446a98SMatthew Knepley 
13883f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
1389e2446a98SMatthew Knepley 
1390e2446a98SMatthew Knepley    Input Parameters:
1391e2446a98SMatthew Knepley +  opt - the option one is seeking
1392e2446a98SMatthew Knepley .  text - short string describing option
1393e2446a98SMatthew Knepley .  man - manual page for option
1394e2446a98SMatthew Knepley -  nmax - maximum number of values
1395e2446a98SMatthew Knepley 
1396e2446a98SMatthew Knepley    Output Parameter:
1397e2446a98SMatthew Knepley +  value - location to copy values
1398e2446a98SMatthew Knepley .  nmax - actual number of values found
1399e2446a98SMatthew Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
1400e2446a98SMatthew Knepley 
1401e2446a98SMatthew Knepley    Level: beginner
1402e2446a98SMatthew Knepley 
1403e2446a98SMatthew Knepley    Notes:
1404e2446a98SMatthew Knepley    The user should pass in an array of doubles
1405e2446a98SMatthew Knepley 
1406e2446a98SMatthew Knepley    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
1407e2446a98SMatthew Knepley 
1408e2446a98SMatthew Knepley    Concepts: options database^array of strings
1409e2446a98SMatthew Knepley 
1410e2446a98SMatthew Knepley .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1411acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
1412e2446a98SMatthew Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1413e2446a98SMatthew Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1414acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1415e2446a98SMatthew Knepley           PetscOptionsList(), PetscOptionsEList()
1416e2446a98SMatthew Knepley @*/
14177087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolArray(const char opt[],const char text[],const char man[],PetscBool  value[],PetscInt *n,PetscBool  *set)
1418e2446a98SMatthew Knepley {
1419e2446a98SMatthew Knepley   PetscErrorCode ierr;
1420e2446a98SMatthew Knepley   PetscInt       i;
14211ae3d29cSBarry Smith   PetscOptions   amsopt;
1422e2446a98SMatthew Knepley 
1423e2446a98SMatthew Knepley   PetscFunctionBegin;
14241ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
1425ace3abfcSBarry Smith     PetscBool  *vals;
14261ae3d29cSBarry Smith 
14271ae3d29cSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_LOGICAL_ARRAY,&amsopt);CHKERRQ(ierr);
1428ace3abfcSBarry Smith     ierr = PetscMalloc((*n)*sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1429ace3abfcSBarry Smith     vals = (PetscBool*)amsopt->data;
14301ae3d29cSBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
14311ae3d29cSBarry Smith     amsopt->arraylength = *n;
14321ae3d29cSBarry Smith   }
1433acfcf0e5SJed Brown   ierr = PetscOptionsGetBoolArray(PetscOptionsObject.prefix,opt,value,n,set);CHKERRQ(ierr);
1434e2446a98SMatthew Knepley   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
1435e2446a98SMatthew Knepley     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%d",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,value[0]);CHKERRQ(ierr);
1436e2446a98SMatthew Knepley     for (i=1; i<*n; i++) {
1437e2446a98SMatthew Knepley       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,",%d",value[i]);CHKERRQ(ierr);
1438e2446a98SMatthew Knepley     }
14392aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
1440e2446a98SMatthew Knepley   }
1441e2446a98SMatthew Knepley   PetscFunctionReturn(0);
1442e2446a98SMatthew Knepley }
1443e2446a98SMatthew Knepley 
14448cc676e6SMatthew G Knepley #undef __FUNCT__
14458cc676e6SMatthew G Knepley #define __FUNCT__ "PetscOptionsViewer"
14468cc676e6SMatthew G Knepley /*@C
14478cc676e6SMatthew G Knepley    PetscOptionsInt - Gets a viewer appropriate for the type indicated by the user
14488cc676e6SMatthew G Knepley 
14498cc676e6SMatthew G Knepley    Logically Collective on the communicator passed in PetscOptionsBegin()
14508cc676e6SMatthew G Knepley 
14518cc676e6SMatthew G Knepley    Input Parameters:
14528cc676e6SMatthew G Knepley +  opt - option name
14538cc676e6SMatthew G Knepley .  text - short string that describes the option
14548cc676e6SMatthew G Knepley -  man - manual page with additional information on option
14558cc676e6SMatthew G Knepley 
14568cc676e6SMatthew G Knepley    Output Parameter:
14578cc676e6SMatthew G Knepley +  viewer - the viewer
14588cc676e6SMatthew G Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
14598cc676e6SMatthew G Knepley 
14608cc676e6SMatthew G Knepley    Level: beginner
14618cc676e6SMatthew G Knepley 
14628cc676e6SMatthew G Knepley    Concepts: options database^has int
14638cc676e6SMatthew G Knepley 
14648cc676e6SMatthew G Knepley    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
14658cc676e6SMatthew G Knepley If no value is provided ascii:stdout is used
14668cc676e6SMatthew G Knepley $       ascii[:[filename][:format]]   defaults to stdout - format can be one of info, info_detailed, or matlab, for example ascii::info prints just the info
14678cc676e6SMatthew G Knepley $                                     about the object to standard out
14688cc676e6SMatthew G Knepley $       binary[:filename]   defaults to binaryoutput
14698cc676e6SMatthew G Knepley $       draw
14708cc676e6SMatthew G Knepley $       socket[:port]    defaults to the standard output port
14718cc676e6SMatthew G Knepley 
1472*cffb1e40SBarry Smith    Use PetscRestoreViewerDestroy() after using the viewer, otherwise a memory leak will occur
14738cc676e6SMatthew G Knepley 
14748cc676e6SMatthew G Knepley .seealso: PetscOptionsGetViewer(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
14758cc676e6SMatthew G Knepley           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
14768cc676e6SMatthew G Knepley           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
14778cc676e6SMatthew G Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
14788cc676e6SMatthew G Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
14798cc676e6SMatthew G Knepley           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
14808cc676e6SMatthew G Knepley           PetscOptionsList(), PetscOptionsEList()
14818cc676e6SMatthew G Knepley @*/
1482*cffb1e40SBarry Smith PetscErrorCode  PetscOptionsViewer(const char opt[],const char text[],const char man[],PetscViewer *viewer,PetscViewerFormat *format,PetscBool  *set)
14838cc676e6SMatthew G Knepley {
14848cc676e6SMatthew G Knepley   PetscErrorCode ierr;
14858cc676e6SMatthew G Knepley   PetscOptions   amsopt;
14868cc676e6SMatthew G Knepley 
14878cc676e6SMatthew G Knepley   PetscFunctionBegin;
14888cc676e6SMatthew G Knepley   if (!PetscOptionsPublishCount) {
14898cc676e6SMatthew G Knepley     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
14908cc676e6SMatthew G Knepley     ierr = PetscMalloc(sizeof(PetscInt),&amsopt->data);CHKERRQ(ierr);
14918cc676e6SMatthew G Knepley     *(const char**)amsopt->data = "";
14928cc676e6SMatthew G Knepley   }
1493*cffb1e40SBarry Smith   ierr = PetscOptionsGetViewer(PetscOptionsObject.comm,PetscOptionsObject.prefix,opt,viewer,format,set);CHKERRQ(ierr);
14948cc676e6SMatthew G Knepley   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
14958cc676e6SMatthew G Knepley     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%s>: %s (%s)\n",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,"",text,ManSection(man));CHKERRQ(ierr);
14968cc676e6SMatthew G Knepley   }
14978cc676e6SMatthew G Knepley   PetscFunctionReturn(0);
14988cc676e6SMatthew G Knepley }
14998cc676e6SMatthew G Knepley 
150053acd3b1SBarry Smith 
150153acd3b1SBarry Smith #undef __FUNCT__
150253acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsHead"
150353acd3b1SBarry Smith /*@C
1504b52f573bSBarry Smith      PetscOptionsHead - Puts a heading before listing any more published options. Used, for example,
150553acd3b1SBarry Smith             in KSPSetFromOptions_GMRES().
150653acd3b1SBarry Smith 
15073f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
150853acd3b1SBarry Smith 
150953acd3b1SBarry Smith    Input Parameter:
151053acd3b1SBarry Smith .   head - the heading text
151153acd3b1SBarry Smith 
151253acd3b1SBarry Smith 
151353acd3b1SBarry Smith    Level: intermediate
151453acd3b1SBarry Smith 
151553acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
151653acd3b1SBarry Smith 
1517b52f573bSBarry Smith           Can be followed by a call to PetscOptionsTail() in the same function.
151853acd3b1SBarry Smith 
151953acd3b1SBarry Smith    Concepts: options database^subheading
152053acd3b1SBarry Smith 
152153acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1522acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
152353acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
152453acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1525acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
152653acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
152753acd3b1SBarry Smith @*/
15287087cfbeSBarry Smith PetscErrorCode  PetscOptionsHead(const char head[])
152953acd3b1SBarry Smith {
153053acd3b1SBarry Smith   PetscErrorCode ierr;
153153acd3b1SBarry Smith 
153253acd3b1SBarry Smith   PetscFunctionBegin;
153361b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
153453acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  %s\n",head);CHKERRQ(ierr);
153553acd3b1SBarry Smith   }
153653acd3b1SBarry Smith   PetscFunctionReturn(0);
153753acd3b1SBarry Smith }
153853acd3b1SBarry Smith 
153953acd3b1SBarry Smith 
154053acd3b1SBarry Smith 
154153acd3b1SBarry Smith 
154253acd3b1SBarry Smith 
154353acd3b1SBarry Smith 
1544