xref: /petsc/src/sys/objects/aoptions.c (revision 0298fd7132830bec7daee99a80be0eddb2b310a5)
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;
37a297a907SKarl Rupp 
38c31cb41cSBarry Smith   ierr = PetscFree(PetscOptionsObject.prefix);CHKERRQ(ierr);
3953acd3b1SBarry Smith   ierr = PetscStrallocpy(prefix,&PetscOptionsObject.prefix);CHKERRQ(ierr);
40c31cb41cSBarry Smith   ierr = PetscFree(PetscOptionsObject.title);CHKERRQ(ierr);
4153acd3b1SBarry Smith   ierr = PetscStrallocpy(title,&PetscOptionsObject.title);CHKERRQ(ierr);
4253acd3b1SBarry Smith 
43*0298fd71SBarry Smith   ierr = PetscOptionsHasName(NULL,"-help",&PetscOptionsObject.printhelp);CHKERRQ(ierr);
4453acd3b1SBarry Smith   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1) {
4561b37b28SSatish Balay     if (!PetscOptionsObject.alreadyprinted) {
4653acd3b1SBarry Smith       ierr = (*PetscHelpPrintf)(comm,"%s -------------------------------------------------\n",title);CHKERRQ(ierr);
4753acd3b1SBarry Smith     }
4861b37b28SSatish Balay   }
4953acd3b1SBarry Smith   PetscFunctionReturn(0);
5053acd3b1SBarry Smith }
5153acd3b1SBarry Smith 
523194b578SJed Brown #undef __FUNCT__
533194b578SJed Brown #define __FUNCT__ "PetscObjectOptionsBegin_Private"
543194b578SJed Brown /*
553194b578SJed Brown     Handles setting up the data structure in a call to PetscObjectOptionsBegin()
563194b578SJed Brown */
573194b578SJed Brown PetscErrorCode PetscObjectOptionsBegin_Private(PetscObject obj)
583194b578SJed Brown {
593194b578SJed Brown   PetscErrorCode ierr;
603194b578SJed Brown   char           title[256];
613194b578SJed Brown   PetscBool      flg;
623194b578SJed Brown 
633194b578SJed Brown   PetscFunctionBegin;
643194b578SJed Brown   PetscValidHeader(obj,1);
653194b578SJed Brown   PetscOptionsObject.object         = obj;
663194b578SJed Brown   PetscOptionsObject.alreadyprinted = obj->optionsprinted;
67a297a907SKarl Rupp 
683194b578SJed Brown   ierr = PetscStrcmp(obj->description,obj->class_name,&flg);CHKERRQ(ierr);
693194b578SJed Brown   if (flg) {
708caf3d72SBarry Smith     ierr = PetscSNPrintf(title,sizeof(title),"%s options",obj->class_name);CHKERRQ(ierr);
713194b578SJed Brown   } else {
728caf3d72SBarry Smith     ierr = PetscSNPrintf(title,sizeof(title),"%s (%s) options",obj->description,obj->class_name);CHKERRQ(ierr);
733194b578SJed Brown   }
743194b578SJed Brown   ierr = PetscOptionsBegin_Private(obj->comm,obj->prefix,title,obj->mansec);CHKERRQ(ierr);
753194b578SJed Brown   PetscFunctionReturn(0);
763194b578SJed Brown }
773194b578SJed Brown 
7853acd3b1SBarry Smith /*
7953acd3b1SBarry Smith      Handles adding another option to the list of options within this particular PetscOptionsBegin() PetscOptionsEnd()
8053acd3b1SBarry Smith */
8153acd3b1SBarry Smith #undef __FUNCT__
8253acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsCreate_Private"
83e3ed6ec8SBarry Smith static int PetscOptionsCreate_Private(const char opt[],const char text[],const char man[],PetscOptionType t,PetscOptions *amsopt)
8453acd3b1SBarry Smith {
8553acd3b1SBarry Smith   int          ierr;
8653acd3b1SBarry Smith   PetscOptions next;
873be6e4c3SJed Brown   PetscBool    valid;
8853acd3b1SBarry Smith 
8953acd3b1SBarry Smith   PetscFunctionBegin;
903be6e4c3SJed Brown   ierr = PetscOptionsValidKey(opt,&valid);CHKERRQ(ierr);
913be6e4c3SJed Brown   if (!valid) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"The option '%s' is not a valid key",opt);
923be6e4c3SJed Brown 
936e655a9bSJed Brown   ierr            = PetscNew(struct _n_PetscOptions,amsopt);CHKERRQ(ierr);
9453acd3b1SBarry Smith   (*amsopt)->next = 0;
9553acd3b1SBarry Smith   (*amsopt)->set  = PETSC_FALSE;
966356e834SBarry Smith   (*amsopt)->type = t;
9753acd3b1SBarry Smith   (*amsopt)->data = 0;
9861b37b28SSatish Balay 
9953acd3b1SBarry Smith   ierr = PetscStrallocpy(text,&(*amsopt)->text);CHKERRQ(ierr);
10053acd3b1SBarry Smith   ierr = PetscStrallocpy(opt,&(*amsopt)->option);CHKERRQ(ierr);
1016356e834SBarry Smith   ierr = PetscStrallocpy(man,&(*amsopt)->man);CHKERRQ(ierr);
10253acd3b1SBarry Smith 
103a297a907SKarl Rupp   if (!PetscOptionsObject.next) PetscOptionsObject.next = *amsopt;
104a297a907SKarl Rupp   else {
10553acd3b1SBarry Smith     next = PetscOptionsObject.next;
10653acd3b1SBarry Smith     while (next->next) next = next->next;
10753acd3b1SBarry Smith     next->next = *amsopt;
10853acd3b1SBarry Smith   }
10953acd3b1SBarry Smith   PetscFunctionReturn(0);
11053acd3b1SBarry Smith }
11153acd3b1SBarry Smith 
11253acd3b1SBarry Smith #undef __FUNCT__
113aee2cecaSBarry Smith #define __FUNCT__ "PetscScanString"
114aee2cecaSBarry Smith /*
1153fc1eb6aSBarry Smith     PetscScanString -  Gets user input via stdin from process and broadcasts to all processes
1163fc1eb6aSBarry Smith 
1173fc1eb6aSBarry Smith     Collective on MPI_Comm
1183fc1eb6aSBarry Smith 
1193fc1eb6aSBarry Smith    Input Parameters:
1203fc1eb6aSBarry Smith +     commm - communicator for the broadcast, must be PETSC_COMM_WORLD
1213fc1eb6aSBarry Smith .     n - length of the string, must be the same on all processes
1223fc1eb6aSBarry Smith -     str - location to store input
123aee2cecaSBarry Smith 
124aee2cecaSBarry Smith     Bugs:
125aee2cecaSBarry Smith .   Assumes process 0 of the given communicator has access to stdin
126aee2cecaSBarry Smith 
127aee2cecaSBarry Smith */
1283fc1eb6aSBarry Smith static PetscErrorCode PetscScanString(MPI_Comm comm,size_t n,char str[])
129aee2cecaSBarry Smith {
130330cf3c9SBarry Smith   size_t         i;
131aee2cecaSBarry Smith   char           c;
1323fc1eb6aSBarry Smith   PetscMPIInt    rank,nm;
133aee2cecaSBarry Smith   PetscErrorCode ierr;
134aee2cecaSBarry Smith 
135aee2cecaSBarry Smith   PetscFunctionBegin;
136aee2cecaSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
137aee2cecaSBarry Smith   if (!rank) {
138aee2cecaSBarry Smith     c = (char) getchar();
139aee2cecaSBarry Smith     i = 0;
140aee2cecaSBarry Smith     while (c != '\n' && i < n-1) {
141aee2cecaSBarry Smith       str[i++] = c;
142aee2cecaSBarry Smith       c = (char)getchar();
143aee2cecaSBarry Smith     }
144aee2cecaSBarry Smith     str[i] = 0;
145aee2cecaSBarry Smith   }
1464dc2109aSBarry Smith   ierr = PetscMPIIntCast(n,&nm);CHKERRQ(ierr);
1473fc1eb6aSBarry Smith   ierr = MPI_Bcast(str,nm,MPI_CHAR,0,comm);CHKERRQ(ierr);
148aee2cecaSBarry Smith   PetscFunctionReturn(0);
149aee2cecaSBarry Smith }
150aee2cecaSBarry Smith 
151aee2cecaSBarry Smith #undef __FUNCT__
152aee2cecaSBarry Smith #define __FUNCT__ "PetscOptionsGetFromTextInput"
153aee2cecaSBarry Smith /*
1543cc1e11dSBarry Smith     PetscOptionsGetFromTextInput - Presents all the PETSc Options processed by the program so the user may change them at runtime
155aee2cecaSBarry Smith 
156aee2cecaSBarry Smith     Notes: this isn't really practical, it is just to demonstrate the principle
157aee2cecaSBarry Smith 
158aee2cecaSBarry Smith     Bugs:
159aee2cecaSBarry Smith +    All processes must traverse through the exact same set of option queries do to the call to PetscScanString()
1603cc1e11dSBarry Smith .    Internal strings have arbitrary length and string copies are not checked that they fit into string space
161aee2cecaSBarry Smith -    Only works for PetscInt == int, PetscReal == double etc
162aee2cecaSBarry Smith 
1633cc1e11dSBarry Smith     Developer Notes: Normally the GUI that presents the options the user and retrieves the values would be running in a different
1643cc1e11dSBarry Smith      address space and communicating with the PETSc program
1653cc1e11dSBarry Smith 
166aee2cecaSBarry Smith */
167aee2cecaSBarry Smith PetscErrorCode PetscOptionsGetFromTextInput()
1686356e834SBarry Smith {
1696356e834SBarry Smith   PetscErrorCode ierr;
1706356e834SBarry Smith   PetscOptions   next = PetscOptionsObject.next;
1716356e834SBarry Smith   char           str[512];
172a4404d99SBarry Smith   PetscInt       id;
173a4404d99SBarry Smith   PetscReal      ir,*valr;
174330cf3c9SBarry Smith   PetscInt       *vald;
175330cf3c9SBarry Smith   size_t         i;
1766356e834SBarry Smith 
177e26ddf31SBarry Smith   ierr = (*PetscPrintf)(PETSC_COMM_WORLD,"%s -------------------------------------------------\n",PetscOptionsObject.title);CHKERRQ(ierr);
1786356e834SBarry Smith   while (next) {
1796356e834SBarry Smith     switch (next->type) {
1806356e834SBarry Smith     case OPTION_HEAD:
1816356e834SBarry Smith       break;
182e26ddf31SBarry Smith     case OPTION_INT_ARRAY:
183e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",next->option+1);CHKERRQ(ierr);
184e26ddf31SBarry Smith       vald = (PetscInt*) next->data;
185e26ddf31SBarry Smith       for (i=0; i<next->arraylength; i++) {
186e26ddf31SBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"%d",vald[i]);CHKERRQ(ierr);
187e26ddf31SBarry Smith         if (i < next->arraylength-1) {
188e26ddf31SBarry Smith           ierr = PetscPrintf(PETSC_COMM_WORLD,",");CHKERRQ(ierr);
189e26ddf31SBarry Smith         }
190e26ddf31SBarry Smith       }
191e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,">: %s (%s)",next->text,next->man);CHKERRQ(ierr);
192e26ddf31SBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
193e26ddf31SBarry Smith       if (str[0]) {
194e26ddf31SBarry Smith         PetscToken token;
195e26ddf31SBarry Smith         PetscInt   n=0,nmax = next->arraylength,*dvalue = (PetscInt*)next->data,start,end;
196e26ddf31SBarry Smith         size_t     len;
197e26ddf31SBarry Smith         char       *value;
198ace3abfcSBarry Smith         PetscBool  foundrange;
199e26ddf31SBarry Smith 
200e26ddf31SBarry Smith         next->set = PETSC_TRUE;
201e26ddf31SBarry Smith         value     = str;
202e26ddf31SBarry Smith         ierr      = PetscTokenCreate(value,',',&token);CHKERRQ(ierr);
203e26ddf31SBarry Smith         ierr      = PetscTokenFind(token,&value);CHKERRQ(ierr);
204e26ddf31SBarry Smith         while (n < nmax) {
205e26ddf31SBarry Smith           if (!value) break;
206e26ddf31SBarry Smith 
207e26ddf31SBarry Smith           /* look for form  d-D where d and D are integers */
208e26ddf31SBarry Smith           foundrange = PETSC_FALSE;
209e26ddf31SBarry Smith           ierr       = PetscStrlen(value,&len);CHKERRQ(ierr);
210e26ddf31SBarry Smith           if (value[0] == '-') i=2;
211e26ddf31SBarry Smith           else i=1;
212330cf3c9SBarry Smith           for (;i<len; i++) {
213e26ddf31SBarry Smith             if (value[i] == '-') {
214e32f2f54SBarry Smith               if (i == len-1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Error in %D-th array entry %s\n",n,value);
215e26ddf31SBarry Smith               value[i] = 0;
216cfbddea1SSatish Balay               ierr     = PetscOptionsStringToInt(value,&start);CHKERRQ(ierr);
217cfbddea1SSatish Balay               ierr     = PetscOptionsStringToInt(value+i+1,&end);CHKERRQ(ierr);
218e32f2f54SBarry 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);
219e32f2f54SBarry 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);
220e26ddf31SBarry Smith               for (; start<end; start++) {
221e26ddf31SBarry Smith                 *dvalue = start; dvalue++;n++;
222e26ddf31SBarry Smith               }
223e26ddf31SBarry Smith               foundrange = PETSC_TRUE;
224e26ddf31SBarry Smith               break;
225e26ddf31SBarry Smith             }
226e26ddf31SBarry Smith           }
227e26ddf31SBarry Smith           if (!foundrange) {
228cfbddea1SSatish Balay             ierr = PetscOptionsStringToInt(value,dvalue);CHKERRQ(ierr);
229e26ddf31SBarry Smith             dvalue++;
230e26ddf31SBarry Smith             n++;
231e26ddf31SBarry Smith           }
232e26ddf31SBarry Smith           ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
233e26ddf31SBarry Smith         }
2348c74ee41SBarry Smith         ierr = PetscTokenDestroy(&token);CHKERRQ(ierr);
235e26ddf31SBarry Smith       }
236e26ddf31SBarry Smith       break;
237e26ddf31SBarry Smith     case OPTION_REAL_ARRAY:
238e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",next->option+1);CHKERRQ(ierr);
239e26ddf31SBarry Smith       valr = (PetscReal*) next->data;
240e26ddf31SBarry Smith       for (i=0; i<next->arraylength; i++) {
241e26ddf31SBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"%g",valr[i]);CHKERRQ(ierr);
242e26ddf31SBarry Smith         if (i < next->arraylength-1) {
243e26ddf31SBarry Smith           ierr = PetscPrintf(PETSC_COMM_WORLD,",");CHKERRQ(ierr);
244e26ddf31SBarry Smith         }
245e26ddf31SBarry Smith       }
246e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,">: %s (%s)",next->text,next->man);CHKERRQ(ierr);
247e26ddf31SBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
248e26ddf31SBarry Smith       if (str[0]) {
249e26ddf31SBarry Smith         PetscToken token;
250e26ddf31SBarry Smith         PetscInt   n = 0,nmax = next->arraylength;
251e26ddf31SBarry Smith         PetscReal  *dvalue = (PetscReal*)next->data;
252e26ddf31SBarry Smith         char       *value;
253e26ddf31SBarry Smith 
254e26ddf31SBarry Smith         next->set = PETSC_TRUE;
255e26ddf31SBarry Smith         value     = str;
256e26ddf31SBarry Smith         ierr      = PetscTokenCreate(value,',',&token);CHKERRQ(ierr);
257e26ddf31SBarry Smith         ierr      = PetscTokenFind(token,&value);CHKERRQ(ierr);
258e26ddf31SBarry Smith         while (n < nmax) {
259e26ddf31SBarry Smith           if (!value) break;
260cfbddea1SSatish Balay           ierr = PetscOptionsStringToReal(value,dvalue);CHKERRQ(ierr);
261e26ddf31SBarry Smith           dvalue++;
262e26ddf31SBarry Smith           n++;
263e26ddf31SBarry Smith           ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
264e26ddf31SBarry Smith         }
2658c74ee41SBarry Smith         ierr = PetscTokenDestroy(&token);CHKERRQ(ierr);
266e26ddf31SBarry Smith       }
267e26ddf31SBarry Smith       break;
2686356e834SBarry Smith     case OPTION_INT:
269e26ddf31SBarry 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);
2703fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
2713fc1eb6aSBarry Smith       if (str[0]) {
272c272547aSJed Brown #if defined(PETSC_USE_64BIT_INDICES)
273c272547aSJed Brown         sscanf(str,"%lld",&id);
274c272547aSJed Brown #else
275aee2cecaSBarry Smith         sscanf(str,"%d",&id);
276c272547aSJed Brown #endif
277aee2cecaSBarry Smith         next->set = PETSC_TRUE;
278a297a907SKarl Rupp 
279aee2cecaSBarry Smith         *((PetscInt*)next->data) = id;
280aee2cecaSBarry Smith       }
281aee2cecaSBarry Smith       break;
282aee2cecaSBarry Smith     case OPTION_REAL:
283e26ddf31SBarry 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);
2843fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
2853fc1eb6aSBarry Smith       if (str[0]) {
286ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
287a4404d99SBarry Smith         sscanf(str,"%e",&ir);
288ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL_DOUBLE)
289aee2cecaSBarry Smith         sscanf(str,"%le",&ir);
290ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL___FLOAT128)
291d9822059SBarry Smith         ir = strtoflt128(str,0);
292d9822059SBarry Smith #else
293513dbe71SLisandro Dalcin         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Unknown scalar type");
294a4404d99SBarry Smith #endif
295aee2cecaSBarry Smith         next->set                 = PETSC_TRUE;
296aee2cecaSBarry Smith         *((PetscReal*)next->data) = ir;
297aee2cecaSBarry Smith       }
298aee2cecaSBarry Smith       break;
299aee2cecaSBarry Smith     case OPTION_LOGICAL:
300aee2cecaSBarry Smith     case OPTION_STRING:
301e26ddf31SBarry 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);
3023fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3033fc1eb6aSBarry Smith       if (str[0]) {
304aee2cecaSBarry Smith         next->set = PETSC_TRUE;
305a297a907SKarl Rupp 
306330cf3c9SBarry Smith         ierr = PetscStrcpy((char*)next->data,str);CHKERRQ(ierr);
3076356e834SBarry Smith       }
3086356e834SBarry Smith       break;
3093cc1e11dSBarry Smith     case OPTION_LIST:
310140e18c1SBarry Smith       ierr = PetscFunctionListPrintTypes(PETSC_COMM_WORLD,stdout,PetscOptionsObject.prefix,next->option,next->text,next->man,next->flist,(char*)next->data);CHKERRQ(ierr);
3113cc1e11dSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3123cc1e11dSBarry Smith       if (str[0]) {
3133cc1e11dSBarry Smith         PetscOptionsObject.changedmethod = PETSC_TRUE;
3143cc1e11dSBarry Smith         next->set = PETSC_TRUE;
315330cf3c9SBarry Smith         ierr = PetscStrcpy((char*)next->data,str);CHKERRQ(ierr);
3163cc1e11dSBarry Smith       }
3173cc1e11dSBarry Smith       break;
318b432afdaSMatthew Knepley     default:
319b432afdaSMatthew Knepley       break;
3206356e834SBarry Smith     }
3216356e834SBarry Smith     next = next->next;
3226356e834SBarry Smith   }
3236356e834SBarry Smith   PetscFunctionReturn(0);
3246356e834SBarry Smith }
3256356e834SBarry Smith 
326b3506946SBarry Smith #if defined(PETSC_HAVE_AMS)
3271ae3d29cSBarry Smith #define CHKERRAMS(err)  if (err) {char *msg; AMS_Explain_error((err), &(msg)); SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"AMS Error: %s",msg);}
3281ae3d29cSBarry 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);}
329d5649816SBarry Smith 
330d5649816SBarry Smith static int count = 0;
331d5649816SBarry Smith 
332b3506946SBarry Smith #undef __FUNCT__
333d5649816SBarry Smith #define __FUNCT__ "PetscOptionsAMSDestroy"
334d5649816SBarry Smith PetscErrorCode PetscOptionsAMSDestroy(void)
335d5649816SBarry Smith {
336d5649816SBarry Smith   PetscErrorCode ierr;
337d5649816SBarry Smith   AMS_Comm       acomm = -1;
338d5649816SBarry Smith   AMS_Memory     amem  = -1;
339d5649816SBarry Smith   char           options[16];
340d5649816SBarry Smith   const char     *string = "Exit";
341d5649816SBarry Smith 
342d5649816SBarry Smith   /* the next line is a bug, this will only work if all processors are here, the comm passed in is ignored!!! */
343d5649816SBarry Smith   ierr = PetscViewerAMSGetAMSComm(PETSC_VIEWER_AMS_(PETSC_COMM_WORLD),&acomm);CHKERRQ(ierr);
344d5649816SBarry Smith   sprintf(options,"Options_%d",count++);
345d5649816SBarry Smith   ierr = AMS_Memory_create(acomm,options,&amem);CHKERRAMS(ierr);
346d5649816SBarry Smith   ierr = AMS_Memory_add_field(amem,"Exit",&string,1,AMS_STRING,AMS_READ,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,"Exit");
347d5649816SBarry Smith 
348d5649816SBarry Smith   ierr = AMS_Memory_take_access(amem);CHKERRAMS(ierr);
349d5649816SBarry Smith   ierr = AMS_Memory_publish(amem);CHKERRAMS(ierr);
350d5649816SBarry Smith   ierr = AMS_Memory_grant_access(amem);CHKERRAMS(ierr);
351d5649816SBarry Smith   /* wait until accessor has unlocked the memory */
352d5649816SBarry Smith   ierr = AMS_Memory_lock(amem,0);CHKERRAMS(ierr);
353d5649816SBarry Smith   ierr = AMS_Memory_take_access(amem);CHKERRAMS(ierr);
354d5649816SBarry Smith   ierr = AMS_Memory_grant_access(amem);CHKERRAMS(ierr);
355d5649816SBarry Smith   ierr = AMS_Memory_destroy(amem);CHKERRAMS(ierr);
356d5649816SBarry Smith   PetscFunctionReturn(0);
357d5649816SBarry Smith }
358d5649816SBarry Smith 
359d5649816SBarry Smith #undef __FUNCT__
360d5649816SBarry Smith #define __FUNCT__ "PetscOptionsAMSInput"
361b3506946SBarry Smith /*
362d5649816SBarry Smith     PetscOptionsAMSInput - Presents all the PETSc Options processed by the program so the user may change them at runtime using the AMS
363b3506946SBarry Smith 
364b3506946SBarry Smith     Bugs:
365b3506946SBarry Smith +    All processes must traverse through the exact same set of option queries do to the call to PetscScanString()
366b3506946SBarry Smith .    Internal strings have arbitrary length and string copies are not checked that they fit into string space
367b3506946SBarry Smith -    Only works for PetscInt == int, PetscReal == double etc
368b3506946SBarry Smith 
369b3506946SBarry Smith 
370b3506946SBarry Smith */
371d5649816SBarry Smith PetscErrorCode PetscOptionsAMSInput()
372b3506946SBarry Smith {
373b3506946SBarry Smith   PetscErrorCode ierr;
374b3506946SBarry Smith   PetscOptions   next     = PetscOptionsObject.next;
375d5649816SBarry Smith   static int     mancount = 0;
376b3506946SBarry Smith   char           options[16];
377b3506946SBarry Smith   AMS_Comm       acomm         = -1;
378b3506946SBarry Smith   AMS_Memory     amem          = -1;
379ace3abfcSBarry Smith   PetscBool      changedmethod = PETSC_FALSE;
3809f32e415SBarry Smith   char           manname[16];
381b3506946SBarry Smith 
382b3506946SBarry Smith   /* the next line is a bug, this will only work if all processors are here, the comm passed in is ignored!!! */
383b3506946SBarry Smith   ierr = PetscViewerAMSGetAMSComm(PETSC_VIEWER_AMS_(PETSC_COMM_WORLD),&acomm);CHKERRQ(ierr);
384b3506946SBarry Smith   sprintf(options,"Options_%d",count++);
3851ae3d29cSBarry Smith   ierr = AMS_Memory_create(acomm,options,&amem);CHKERRAMS(ierr);
3861ae3d29cSBarry Smith   ierr = AMS_Memory_take_access(amem);CHKERRAMS(ierr);
387a297a907SKarl Rupp 
3881bc75a8dSBarry Smith   PetscOptionsObject.pprefix = PetscOptionsObject.prefix; /* AMS will change this, so cannot pass prefix directly */
3891bc75a8dSBarry Smith 
3901ae3d29cSBarry 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);
3910c74a584SJed Brown   /* ierr = AMS_Memory_add_field(amem,"mansec",&PetscOptionsObject.pprefix,1,AMS_STRING,AMS_READ,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMS(ierr); */
3921ae3d29cSBarry Smith   ierr = AMS_Memory_add_field(amem,"ChangedMethod",&changedmethod,1,AMS_BOOLEAN,AMS_WRITE,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,"ChangedMethod");
393b3506946SBarry Smith 
394b3506946SBarry Smith   while (next) {
3951ae3d29cSBarry 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);
3961bc75a8dSBarry Smith     ierr =  PetscMalloc(sizeof(char*),&next->pman);CHKERRQ(ierr);
397a297a907SKarl Rupp 
3981bc75a8dSBarry Smith     *(char**)next->pman = next->man;
3999f32e415SBarry Smith     sprintf(manname,"man_%d",mancount++);
4001ae3d29cSBarry Smith     ierr = AMS_Memory_add_field(amem,manname,next->pman,1,AMS_STRING,AMS_READ,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,manname);
4019f32e415SBarry Smith 
402b3506946SBarry Smith     switch (next->type) {
403b3506946SBarry Smith     case OPTION_HEAD:
404b3506946SBarry Smith       break;
405b3506946SBarry Smith     case OPTION_INT_ARRAY:
4061ae3d29cSBarry 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);
407b3506946SBarry Smith       break;
408b3506946SBarry Smith     case OPTION_REAL_ARRAY:
4091ae3d29cSBarry 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);
410b3506946SBarry Smith       break;
411b3506946SBarry Smith     case OPTION_INT:
4121ae3d29cSBarry 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);
413b3506946SBarry Smith       break;
414b3506946SBarry Smith     case OPTION_REAL:
4151ae3d29cSBarry 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);
416b3506946SBarry Smith       break;
417b3506946SBarry Smith     case OPTION_LOGICAL:
4181ae3d29cSBarry 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);
4191ae3d29cSBarry Smith       break;
4201ae3d29cSBarry Smith     case OPTION_LOGICAL_ARRAY:
4211ae3d29cSBarry 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);
42271f08665SBarry Smith       break;
423b3506946SBarry Smith     case OPTION_STRING:
4241ae3d29cSBarry 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);
4251ae3d29cSBarry Smith       break;
4261ae3d29cSBarry Smith     case OPTION_STRING_ARRAY:
4271ae3d29cSBarry 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);
428b3506946SBarry Smith       break;
429b3506946SBarry Smith     case OPTION_LIST:
43071f08665SBarry Smith       {PetscInt ntext;
43171f08665SBarry Smith       char      ldefault[128];
43271f08665SBarry Smith       ierr = PetscStrcpy(ldefault,"DEFAULT:");CHKERRQ(ierr);
43371f08665SBarry Smith       ierr = PetscStrcat(ldefault,next->text);CHKERRQ(ierr);
4341ae3d29cSBarry Smith       ierr = AMS_Memory_add_field(amem,ldefault,next->data,1,AMS_STRING,AMS_WRITE,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,ldefault);
435140e18c1SBarry Smith       ierr = PetscFunctionListGet(next->flist,(const char***)&next->edata,&ntext);CHKERRQ(ierr);
436d5649816SBarry 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);
4371ae3d29cSBarry Smith       break;}
4381ae3d29cSBarry Smith     case OPTION_ELIST:
439d5649816SBarry Smith       {PetscInt ntext = next->nlist;
4401ae3d29cSBarry Smith       char      ldefault[128];
4411ae3d29cSBarry Smith       ierr = PetscStrcpy(ldefault,"DEFAULT:");CHKERRQ(ierr);
4421ae3d29cSBarry Smith       ierr = PetscStrcat(ldefault,next->text);CHKERRQ(ierr);
4431ae3d29cSBarry Smith       ierr = AMS_Memory_add_field(amem,ldefault,next->data,1,AMS_STRING,AMS_WRITE,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,ldefault);
4441ae3d29cSBarry Smith       ierr = PetscMalloc((ntext+1)*sizeof(char**),&next->edata);CHKERRQ(ierr);
4451ae3d29cSBarry Smith       ierr = PetscMemcpy(next->edata,next->list,ntext*sizeof(char*));CHKERRQ(ierr);
4461ae3d29cSBarry 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);
44771f08665SBarry Smith       break;}
448b3506946SBarry Smith     default:
449b3506946SBarry Smith       break;
450b3506946SBarry Smith     }
451b3506946SBarry Smith     next = next->next;
452b3506946SBarry Smith   }
453b3506946SBarry Smith 
4541ae3d29cSBarry Smith   ierr = AMS_Memory_publish(amem);CHKERRAMS(ierr);
4551ae3d29cSBarry Smith   ierr = AMS_Memory_grant_access(amem);CHKERRAMS(ierr);
456b3506946SBarry Smith   /* wait until accessor has unlocked the memory */
4571ae3d29cSBarry Smith   ierr = AMS_Memory_lock(amem,0);CHKERRAMS(ierr);
4581ae3d29cSBarry Smith   ierr = AMS_Memory_take_access(amem);CHKERRAMS(ierr);
459b3506946SBarry Smith 
460b3506946SBarry Smith   /* reset counter to -2; this updates the screen with the new options for the selected method */
461b3506946SBarry Smith   if (changedmethod) PetscOptionsPublishCount = -2;
462b3506946SBarry Smith 
4631ae3d29cSBarry Smith   ierr = AMS_Memory_grant_access(amem);CHKERRAMS(ierr);
4641ae3d29cSBarry Smith   ierr = AMS_Memory_destroy(amem);CHKERRAMS(ierr);
465b3506946SBarry Smith   PetscFunctionReturn(0);
466b3506946SBarry Smith }
467b3506946SBarry Smith #endif
468b3506946SBarry Smith 
4696356e834SBarry Smith #undef __FUNCT__
47053acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEnd_Private"
47153acd3b1SBarry Smith PetscErrorCode PetscOptionsEnd_Private(void)
47253acd3b1SBarry Smith {
47353acd3b1SBarry Smith   PetscErrorCode ierr;
4746356e834SBarry Smith   PetscOptions   last;
4756356e834SBarry Smith   char           option[256],value[1024],tmp[32];
476330cf3c9SBarry Smith   size_t         j;
47753acd3b1SBarry Smith 
47853acd3b1SBarry Smith   PetscFunctionBegin;
4791bc75a8dSBarry Smith   CHKMEMQ;
480aee2cecaSBarry Smith   if (PetscOptionsObject.next) {
481b3506946SBarry Smith     if (!PetscOptionsPublishCount) {
482b3506946SBarry Smith #if defined(PETSC_HAVE_AMS)
483d5649816SBarry Smith       ierr = PetscOptionsAMSInput();CHKERRQ(ierr);
484b3506946SBarry Smith #else
48571f08665SBarry Smith       ierr = PetscOptionsGetFromTextInput();CHKERRQ(ierr);
486b3506946SBarry Smith #endif
487aee2cecaSBarry Smith     }
488aee2cecaSBarry Smith   }
4896356e834SBarry Smith 
490c31cb41cSBarry Smith   ierr = PetscFree(PetscOptionsObject.title);CHKERRQ(ierr);
491c31cb41cSBarry Smith   ierr = PetscFree(PetscOptionsObject.prefix);CHKERRQ(ierr);
4926356e834SBarry Smith 
4936356e834SBarry Smith   /* reset counter to -2; this updates the screen with the new options for the selected method */
4946356e834SBarry Smith   if (PetscOptionsObject.changedmethod) PetscOptionsPublishCount = -2;
49561b37b28SSatish Balay   /* reset alreadyprinted flag */
49661b37b28SSatish Balay   PetscOptionsObject.alreadyprinted = PETSC_FALSE;
4973194b578SJed Brown   if (PetscOptionsObject.object) PetscOptionsObject.object->optionsprinted = PETSC_TRUE;
498*0298fd71SBarry Smith   PetscOptionsObject.object = NULL;
4996356e834SBarry Smith 
5006356e834SBarry Smith   while (PetscOptionsObject.next) {
5016356e834SBarry Smith     if (PetscOptionsObject.next->set) {
5026356e834SBarry Smith       if (PetscOptionsObject.prefix) {
5036356e834SBarry Smith         ierr = PetscStrcpy(option,"-");CHKERRQ(ierr);
5046356e834SBarry Smith         ierr = PetscStrcat(option,PetscOptionsObject.prefix);CHKERRQ(ierr);
5056356e834SBarry Smith         ierr = PetscStrcat(option,PetscOptionsObject.next->option+1);CHKERRQ(ierr);
5066356e834SBarry Smith       } else {
5076356e834SBarry Smith         ierr = PetscStrcpy(option,PetscOptionsObject.next->option);CHKERRQ(ierr);
5086356e834SBarry Smith       }
5096356e834SBarry Smith 
5106356e834SBarry Smith       switch (PetscOptionsObject.next->type) {
5116356e834SBarry Smith       case OPTION_HEAD:
5126356e834SBarry Smith         break;
513e26ddf31SBarry Smith       case OPTION_INT_ARRAY:
514e26ddf31SBarry Smith         sprintf(value,"%d",(int)((PetscInt*)PetscOptionsObject.next->data)[0]);
515e26ddf31SBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
516e26ddf31SBarry Smith           sprintf(tmp,"%d",(int)((PetscInt*)PetscOptionsObject.next->data)[j]);
517e26ddf31SBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
518e26ddf31SBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
519e26ddf31SBarry Smith         }
520e26ddf31SBarry Smith         break;
5216356e834SBarry Smith       case OPTION_INT:
5227a72a596SBarry Smith         sprintf(value,"%d",(int) *(PetscInt*)PetscOptionsObject.next->data);
5236356e834SBarry Smith         break;
5246356e834SBarry Smith       case OPTION_REAL:
5257a72a596SBarry Smith         sprintf(value,"%g",(double) *(PetscReal*)PetscOptionsObject.next->data);
5266356e834SBarry Smith         break;
5276356e834SBarry Smith       case OPTION_REAL_ARRAY:
5287a72a596SBarry Smith         sprintf(value,"%g",(double)((PetscReal*)PetscOptionsObject.next->data)[0]);
5296356e834SBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
5307a72a596SBarry Smith           sprintf(tmp,"%g",(double)((PetscReal*)PetscOptionsObject.next->data)[j]);
5316356e834SBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5326356e834SBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5336356e834SBarry Smith         }
5346356e834SBarry Smith         break;
5356356e834SBarry Smith       case OPTION_LOGICAL:
53671f08665SBarry Smith         sprintf(value,"%d",*(int*)PetscOptionsObject.next->data);
5376356e834SBarry Smith         break;
5381ae3d29cSBarry Smith       case OPTION_LOGICAL_ARRAY:
539ace3abfcSBarry Smith         sprintf(value,"%d",(int)((PetscBool*)PetscOptionsObject.next->data)[0]);
5401ae3d29cSBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
541ace3abfcSBarry Smith           sprintf(tmp,"%d",(int)((PetscBool*)PetscOptionsObject.next->data)[j]);
5421ae3d29cSBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5431ae3d29cSBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5441ae3d29cSBarry Smith         }
5451ae3d29cSBarry Smith         break;
5466356e834SBarry Smith       case OPTION_LIST:
5471ae3d29cSBarry Smith       case OPTION_ELIST:
54871f08665SBarry Smith         ierr = PetscStrcpy(value,*(char**)PetscOptionsObject.next->data);CHKERRQ(ierr);
5496356e834SBarry Smith         break;
5501ae3d29cSBarry Smith       case OPTION_STRING:
55171f08665SBarry Smith         ierr = PetscStrcpy(value,*(char**)PetscOptionsObject.next->data);CHKERRQ(ierr);
5521ae3d29cSBarry Smith       case OPTION_STRING_ARRAY:
5531ae3d29cSBarry Smith         sprintf(value,"%s",((char**)PetscOptionsObject.next->data)[0]);
5541ae3d29cSBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
5551ae3d29cSBarry Smith           sprintf(tmp,"%s",((char**)PetscOptionsObject.next->data)[j]);
5561ae3d29cSBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5571ae3d29cSBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5581ae3d29cSBarry Smith         }
5596356e834SBarry Smith         break;
5606356e834SBarry Smith       }
5616356e834SBarry Smith       ierr = PetscOptionsSetValue(option,value);CHKERRQ(ierr);
5626356e834SBarry Smith     }
563503cfb0cSBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->text);CHKERRQ(ierr);
564503cfb0cSBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->option);CHKERRQ(ierr);
5656356e834SBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->man);CHKERRQ(ierr);
56605b42c5fSBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->data);CHKERRQ(ierr);
56771f08665SBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->edata);CHKERRQ(ierr);
568a297a907SKarl Rupp 
5696356e834SBarry Smith     last                    = PetscOptionsObject.next;
5706356e834SBarry Smith     PetscOptionsObject.next = PetscOptionsObject.next->next;
5716356e834SBarry Smith     ierr                    = PetscFree(last);CHKERRQ(ierr);
5721bc75a8dSBarry Smith     CHKMEMQ;
5736356e834SBarry Smith   }
5741bc75a8dSBarry Smith   CHKMEMQ;
5756356e834SBarry Smith   PetscOptionsObject.next = 0;
57653acd3b1SBarry Smith   PetscFunctionReturn(0);
57753acd3b1SBarry Smith }
57853acd3b1SBarry Smith 
57953acd3b1SBarry Smith #undef __FUNCT__
58053acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEnum"
58153acd3b1SBarry Smith /*@C
58253acd3b1SBarry Smith    PetscOptionsEnum - Gets the enum value for a particular option in the database.
58353acd3b1SBarry Smith 
5843f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
58553acd3b1SBarry Smith 
58653acd3b1SBarry Smith    Input Parameters:
58753acd3b1SBarry Smith +  opt - option name
58853acd3b1SBarry Smith .  text - short string that describes the option
58953acd3b1SBarry Smith .  man - manual page with additional information on option
59053acd3b1SBarry Smith .  list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null
59153acd3b1SBarry Smith -  defaultv - the default (current) value
59253acd3b1SBarry Smith 
59353acd3b1SBarry Smith    Output Parameter:
59453acd3b1SBarry Smith +  value - the  value to return
59553acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
59653acd3b1SBarry Smith 
59753acd3b1SBarry Smith    Level: beginner
59853acd3b1SBarry Smith 
59953acd3b1SBarry Smith    Concepts: options database
60053acd3b1SBarry Smith 
60153acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
60253acd3b1SBarry Smith 
60353acd3b1SBarry Smith           list is usually something like PCASMTypes or some other predefined list of enum names
60453acd3b1SBarry Smith 
60553acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
606acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
607acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
60853acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
60953acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
610acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
61153acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
61253acd3b1SBarry Smith @*/
6137087cfbeSBarry Smith PetscErrorCode  PetscOptionsEnum(const char opt[],const char text[],const char man[],const char *const *list,PetscEnum defaultv,PetscEnum *value,PetscBool  *set)
61453acd3b1SBarry Smith {
61553acd3b1SBarry Smith   PetscErrorCode ierr;
61653acd3b1SBarry Smith   PetscInt       ntext = 0;
617aa5bb8c0SSatish Balay   PetscInt       tval;
618ace3abfcSBarry Smith   PetscBool      tflg;
61953acd3b1SBarry Smith 
62053acd3b1SBarry Smith   PetscFunctionBegin;
62153acd3b1SBarry Smith   while (list[ntext++]) {
622e32f2f54SBarry Smith     if (ntext > 50) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument appears to be wrong or have more than 50 entries");
62353acd3b1SBarry Smith   }
624e32f2f54SBarry Smith   if (ntext < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument must have at least two entries: typename and type prefix");
62553acd3b1SBarry Smith   ntext -= 3;
626aa5bb8c0SSatish Balay   ierr   = PetscOptionsEList(opt,text,man,list,ntext,list[defaultv],&tval,&tflg);CHKERRQ(ierr);
627aa5bb8c0SSatish Balay   /* with PETSC_USE_64BIT_INDICES sizeof(PetscInt) != sizeof(PetscEnum) */
628aa5bb8c0SSatish Balay   if (tflg) *value = (PetscEnum)tval;
629aa5bb8c0SSatish Balay   if (set)  *set   = tflg;
63053acd3b1SBarry Smith   PetscFunctionReturn(0);
63153acd3b1SBarry Smith }
63253acd3b1SBarry Smith 
63353acd3b1SBarry Smith /* -------------------------------------------------------------------------------------------------------------*/
63453acd3b1SBarry Smith #undef __FUNCT__
63553acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsInt"
63653acd3b1SBarry Smith /*@C
63753acd3b1SBarry Smith    PetscOptionsInt - Gets the integer value for a particular option in the database.
63853acd3b1SBarry Smith 
6393f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
64053acd3b1SBarry Smith 
64153acd3b1SBarry Smith    Input Parameters:
64253acd3b1SBarry Smith +  opt - option name
64353acd3b1SBarry Smith .  text - short string that describes the option
64453acd3b1SBarry Smith .  man - manual page with additional information on option
64553acd3b1SBarry Smith -  defaultv - the default (current) value
64653acd3b1SBarry Smith 
64753acd3b1SBarry Smith    Output Parameter:
64853acd3b1SBarry Smith +  value - the integer value to return
64953acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
65053acd3b1SBarry Smith 
65153acd3b1SBarry Smith    Level: beginner
65253acd3b1SBarry Smith 
65353acd3b1SBarry Smith    Concepts: options database^has int
65453acd3b1SBarry Smith 
65553acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
65653acd3b1SBarry Smith 
65753acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
658acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
659acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
66053acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
66153acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
662acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
66353acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
66453acd3b1SBarry Smith @*/
6657087cfbeSBarry Smith PetscErrorCode  PetscOptionsInt(const char opt[],const char text[],const char man[],PetscInt defaultv,PetscInt *value,PetscBool  *set)
66653acd3b1SBarry Smith {
66753acd3b1SBarry Smith   PetscErrorCode ierr;
6686356e834SBarry Smith   PetscOptions   amsopt;
66953acd3b1SBarry Smith 
67053acd3b1SBarry Smith   PetscFunctionBegin;
671b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
6726356e834SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_INT,&amsopt);CHKERRQ(ierr);
6736356e834SBarry Smith     ierr = PetscMalloc(sizeof(PetscInt),&amsopt->data);CHKERRQ(ierr);
674a297a907SKarl Rupp 
6756356e834SBarry Smith     *(PetscInt*)amsopt->data = defaultv;
676af6d86caSBarry Smith   }
67753acd3b1SBarry Smith   ierr = PetscOptionsGetInt(PetscOptionsObject.prefix,opt,value,set);CHKERRQ(ierr);
67861b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
6792aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%d>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,defaultv,text,ManSection(man));CHKERRQ(ierr);
68053acd3b1SBarry Smith   }
68153acd3b1SBarry Smith   PetscFunctionReturn(0);
68253acd3b1SBarry Smith }
68353acd3b1SBarry Smith 
68453acd3b1SBarry Smith #undef __FUNCT__
68553acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsString"
68653acd3b1SBarry Smith /*@C
68753acd3b1SBarry Smith    PetscOptionsString - Gets the string value for a particular option in the database.
68853acd3b1SBarry Smith 
6893f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
69053acd3b1SBarry Smith 
69153acd3b1SBarry Smith    Input Parameters:
69253acd3b1SBarry Smith +  opt - option name
69353acd3b1SBarry Smith .  text - short string that describes the option
69453acd3b1SBarry Smith .  man - manual page with additional information on option
695bcbf2dc5SJed Brown .  defaultv - the default (current) value
696bcbf2dc5SJed Brown -  len - length of the result string including null terminator
69753acd3b1SBarry Smith 
69853acd3b1SBarry Smith    Output Parameter:
69953acd3b1SBarry Smith +  value - the value to return
70053acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
70153acd3b1SBarry Smith 
70253acd3b1SBarry Smith    Level: beginner
70353acd3b1SBarry Smith 
70453acd3b1SBarry Smith    Concepts: options database^has int
70553acd3b1SBarry Smith 
70653acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
70753acd3b1SBarry Smith 
7087fccdfe4SBarry 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).
7097fccdfe4SBarry Smith 
71053acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
711acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
712acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsReal(), PetscOptionsBool(),
71353acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
71453acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
715acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
71653acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
71753acd3b1SBarry Smith @*/
7187087cfbeSBarry Smith PetscErrorCode  PetscOptionsString(const char opt[],const char text[],const char man[],const char defaultv[],char value[],size_t len,PetscBool  *set)
71953acd3b1SBarry Smith {
72053acd3b1SBarry Smith   PetscErrorCode ierr;
721aee2cecaSBarry Smith   PetscOptions   amsopt;
72253acd3b1SBarry Smith 
72353acd3b1SBarry Smith   PetscFunctionBegin;
724b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
725aee2cecaSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
72671f08665SBarry Smith     ierr = PetscMalloc(sizeof(char*),&amsopt->data);CHKERRQ(ierr);
727a297a907SKarl Rupp 
72871f08665SBarry Smith     *(const char**)amsopt->data = defaultv;
729af6d86caSBarry Smith   }
73053acd3b1SBarry Smith   ierr = PetscOptionsGetString(PetscOptionsObject.prefix,opt,value,len,set);CHKERRQ(ierr);
73161b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
7322aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%s>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,defaultv,text,ManSection(man));CHKERRQ(ierr);
73353acd3b1SBarry Smith   }
73453acd3b1SBarry Smith   PetscFunctionReturn(0);
73553acd3b1SBarry Smith }
73653acd3b1SBarry Smith 
73753acd3b1SBarry Smith #undef __FUNCT__
73853acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsReal"
73953acd3b1SBarry Smith /*@C
74053acd3b1SBarry Smith    PetscOptionsReal - Gets the PetscReal value for a particular option in the database.
74153acd3b1SBarry Smith 
7423f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
74353acd3b1SBarry Smith 
74453acd3b1SBarry Smith    Input Parameters:
74553acd3b1SBarry Smith +  opt - option name
74653acd3b1SBarry Smith .  text - short string that describes the option
74753acd3b1SBarry Smith .  man - manual page with additional information on option
74853acd3b1SBarry Smith -  defaultv - the default (current) value
74953acd3b1SBarry Smith 
75053acd3b1SBarry Smith    Output Parameter:
75153acd3b1SBarry Smith +  value - the value to return
75253acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
75353acd3b1SBarry Smith 
75453acd3b1SBarry Smith    Level: beginner
75553acd3b1SBarry Smith 
75653acd3b1SBarry Smith    Concepts: options database^has int
75753acd3b1SBarry Smith 
75853acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
75953acd3b1SBarry Smith 
76053acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
761acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
762acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
76353acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
76453acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
765acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
76653acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
76753acd3b1SBarry Smith @*/
7687087cfbeSBarry Smith PetscErrorCode  PetscOptionsReal(const char opt[],const char text[],const char man[],PetscReal defaultv,PetscReal *value,PetscBool  *set)
76953acd3b1SBarry Smith {
77053acd3b1SBarry Smith   PetscErrorCode ierr;
771538aa990SBarry Smith   PetscOptions   amsopt;
77253acd3b1SBarry Smith 
77353acd3b1SBarry Smith   PetscFunctionBegin;
774b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
775538aa990SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_REAL,&amsopt);CHKERRQ(ierr);
776538aa990SBarry Smith     ierr = PetscMalloc(sizeof(PetscReal),&amsopt->data);CHKERRQ(ierr);
777a297a907SKarl Rupp 
778538aa990SBarry Smith     *(PetscReal*)amsopt->data = defaultv;
779538aa990SBarry Smith   }
78053acd3b1SBarry Smith   ierr = PetscOptionsGetReal(PetscOptionsObject.prefix,opt,value,set);CHKERRQ(ierr);
78161b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
7822aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%G>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,defaultv,text,ManSection(man));CHKERRQ(ierr);
78353acd3b1SBarry Smith   }
78453acd3b1SBarry Smith   PetscFunctionReturn(0);
78553acd3b1SBarry Smith }
78653acd3b1SBarry Smith 
78753acd3b1SBarry Smith #undef __FUNCT__
78853acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsScalar"
78953acd3b1SBarry Smith /*@C
79053acd3b1SBarry Smith    PetscOptionsScalar - Gets the scalar value for a particular option in the database.
79153acd3b1SBarry Smith 
7923f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
79353acd3b1SBarry Smith 
79453acd3b1SBarry Smith    Input Parameters:
79553acd3b1SBarry Smith +  opt - option name
79653acd3b1SBarry Smith .  text - short string that describes the option
79753acd3b1SBarry Smith .  man - manual page with additional information on option
79853acd3b1SBarry Smith -  defaultv - the default (current) value
79953acd3b1SBarry Smith 
80053acd3b1SBarry Smith    Output Parameter:
80153acd3b1SBarry Smith +  value - the value to return
80253acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
80353acd3b1SBarry Smith 
80453acd3b1SBarry Smith    Level: beginner
80553acd3b1SBarry Smith 
80653acd3b1SBarry Smith    Concepts: options database^has int
80753acd3b1SBarry Smith 
80853acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
80953acd3b1SBarry Smith 
81053acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
811acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
812acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
81353acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
81453acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
815acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
81653acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
81753acd3b1SBarry Smith @*/
8187087cfbeSBarry Smith PetscErrorCode  PetscOptionsScalar(const char opt[],const char text[],const char man[],PetscScalar defaultv,PetscScalar *value,PetscBool  *set)
81953acd3b1SBarry Smith {
82053acd3b1SBarry Smith   PetscErrorCode ierr;
82153acd3b1SBarry Smith 
82253acd3b1SBarry Smith   PetscFunctionBegin;
82353acd3b1SBarry Smith #if !defined(PETSC_USE_COMPLEX)
82453acd3b1SBarry Smith   ierr = PetscOptionsReal(opt,text,man,defaultv,value,set);CHKERRQ(ierr);
82553acd3b1SBarry Smith #else
82653acd3b1SBarry Smith   ierr = PetscOptionsGetScalar(PetscOptionsObject.prefix,opt,value,set);CHKERRQ(ierr);
82753acd3b1SBarry Smith #endif
82853acd3b1SBarry Smith   PetscFunctionReturn(0);
82953acd3b1SBarry Smith }
83053acd3b1SBarry Smith 
83153acd3b1SBarry Smith #undef __FUNCT__
83253acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsName"
83353acd3b1SBarry Smith /*@C
83490d69ab7SBarry 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
83590d69ab7SBarry Smith                       its value is set to false.
83653acd3b1SBarry Smith 
8373f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
83853acd3b1SBarry Smith 
83953acd3b1SBarry Smith    Input Parameters:
84053acd3b1SBarry Smith +  opt - option name
84153acd3b1SBarry Smith .  text - short string that describes the option
84253acd3b1SBarry Smith -  man - manual page with additional information on option
84353acd3b1SBarry Smith 
84453acd3b1SBarry Smith    Output Parameter:
84553acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
84653acd3b1SBarry Smith 
84753acd3b1SBarry Smith    Level: beginner
84853acd3b1SBarry Smith 
84953acd3b1SBarry Smith    Concepts: options database^has int
85053acd3b1SBarry Smith 
85153acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
85253acd3b1SBarry Smith 
85353acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
854acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
855acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
85653acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
85753acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
858acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
85953acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
86053acd3b1SBarry Smith @*/
8617087cfbeSBarry Smith PetscErrorCode  PetscOptionsName(const char opt[],const char text[],const char man[],PetscBool  *flg)
86253acd3b1SBarry Smith {
86353acd3b1SBarry Smith   PetscErrorCode ierr;
8641ae3d29cSBarry Smith   PetscOptions   amsopt;
86553acd3b1SBarry Smith 
86653acd3b1SBarry Smith   PetscFunctionBegin;
8671ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
8681ae3d29cSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_LOGICAL,&amsopt);CHKERRQ(ierr);
869ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
870a297a907SKarl Rupp 
871ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
8721ae3d29cSBarry Smith   }
87353acd3b1SBarry Smith   ierr = PetscOptionsHasName(PetscOptionsObject.prefix,opt,flg);CHKERRQ(ierr);
87461b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
8752aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
87653acd3b1SBarry Smith   }
87753acd3b1SBarry Smith   PetscFunctionReturn(0);
87853acd3b1SBarry Smith }
87953acd3b1SBarry Smith 
88053acd3b1SBarry Smith #undef __FUNCT__
88153acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsList"
88253acd3b1SBarry Smith /*@C
88353acd3b1SBarry Smith      PetscOptionsList - Puts a list of option values that a single one may be selected from
88453acd3b1SBarry Smith 
8853f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
88653acd3b1SBarry Smith 
88753acd3b1SBarry Smith    Input Parameters:
88853acd3b1SBarry Smith +  opt - option name
88953acd3b1SBarry Smith .  text - short string that describes the option
89053acd3b1SBarry Smith .  man - manual page with additional information on option
89153acd3b1SBarry Smith .  list - the possible choices
8923cc1e11dSBarry Smith .  defaultv - the default (current) value
8933cc1e11dSBarry Smith -  len - the length of the character array value
89453acd3b1SBarry Smith 
89553acd3b1SBarry Smith    Output Parameter:
89653acd3b1SBarry Smith +  value - the value to return
89753acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
89853acd3b1SBarry Smith 
89953acd3b1SBarry Smith    Level: intermediate
90053acd3b1SBarry Smith 
90153acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
90253acd3b1SBarry Smith 
90353acd3b1SBarry Smith    See PetscOptionsEList() for when the choices are given in a string array
90453acd3b1SBarry Smith 
90553acd3b1SBarry Smith    To get a listing of all currently specified options,
90688c29154SBarry Smith     see PetscOptionsView() or PetscOptionsGetAll()
90753acd3b1SBarry Smith 
90853acd3b1SBarry Smith    Concepts: options database^list
90953acd3b1SBarry Smith 
91053acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
911acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
91253acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
91353acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
914acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
915171400e9SBarry Smith           PetscOptionsList(), PetscOptionsEList(), PetscOptionsEnum()
91653acd3b1SBarry Smith @*/
917140e18c1SBarry Smith PetscErrorCode  PetscOptionsList(const char opt[],const char ltext[],const char man[],PetscFunctionList list,const char defaultv[],char value[],size_t len,PetscBool  *set)
91853acd3b1SBarry Smith {
91953acd3b1SBarry Smith   PetscErrorCode ierr;
9203cc1e11dSBarry Smith   PetscOptions   amsopt;
92153acd3b1SBarry Smith 
92253acd3b1SBarry Smith   PetscFunctionBegin;
923b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
9243cc1e11dSBarry Smith     ierr = PetscOptionsCreate_Private(opt,ltext,man,OPTION_LIST,&amsopt);CHKERRQ(ierr);
92571f08665SBarry Smith     ierr = PetscMalloc(sizeof(char*),&amsopt->data);CHKERRQ(ierr);
926a297a907SKarl Rupp 
92771f08665SBarry Smith     *(const char**)amsopt->data = defaultv;
928a297a907SKarl Rupp 
9293cc1e11dSBarry Smith     amsopt->flist = list;
9303cc1e11dSBarry Smith   }
93153acd3b1SBarry Smith   ierr = PetscOptionsGetString(PetscOptionsObject.prefix,opt,value,len,set);CHKERRQ(ierr);
93261b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
933140e18c1SBarry Smith     ierr = PetscFunctionListPrintTypes(PetscOptionsObject.comm,stdout,PetscOptionsObject.prefix,opt,ltext,man,list,defaultv);CHKERRQ(ierr);CHKERRQ(ierr);
93453acd3b1SBarry Smith   }
93553acd3b1SBarry Smith   PetscFunctionReturn(0);
93653acd3b1SBarry Smith }
93753acd3b1SBarry Smith 
93853acd3b1SBarry Smith #undef __FUNCT__
93953acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEList"
94053acd3b1SBarry Smith /*@C
94153acd3b1SBarry Smith      PetscOptionsEList - Puts a list of option values that a single one may be selected from
94253acd3b1SBarry Smith 
9433f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
94453acd3b1SBarry Smith 
94553acd3b1SBarry Smith    Input Parameters:
94653acd3b1SBarry Smith +  opt - option name
94753acd3b1SBarry Smith .  ltext - short string that describes the option
94853acd3b1SBarry Smith .  man - manual page with additional information on option
94953acd3b1SBarry Smith .  list - the possible choices
95053acd3b1SBarry Smith .  ntext - number of choices
95153acd3b1SBarry Smith -  defaultv - the default (current) value
95253acd3b1SBarry Smith 
95353acd3b1SBarry Smith    Output Parameter:
95453acd3b1SBarry Smith +  value - the index of the value to return
95553acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
95653acd3b1SBarry Smith 
95753acd3b1SBarry Smith    Level: intermediate
95853acd3b1SBarry Smith 
95953acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
96053acd3b1SBarry Smith 
961140e18c1SBarry Smith    See PetscOptionsList() for when the choices are given in a PetscFunctionList()
96253acd3b1SBarry Smith 
96353acd3b1SBarry Smith    Concepts: options database^list
96453acd3b1SBarry Smith 
96553acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
966acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
96753acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
96853acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
969acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
970171400e9SBarry Smith           PetscOptionsList(), PetscOptionsEnum()
97153acd3b1SBarry Smith @*/
9727087cfbeSBarry 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)
97353acd3b1SBarry Smith {
97453acd3b1SBarry Smith   PetscErrorCode ierr;
97553acd3b1SBarry Smith   PetscInt       i;
9761ae3d29cSBarry Smith   PetscOptions   amsopt;
97753acd3b1SBarry Smith 
97853acd3b1SBarry Smith   PetscFunctionBegin;
9791ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
980d5649816SBarry Smith     ierr = PetscOptionsCreate_Private(opt,ltext,man,OPTION_ELIST,&amsopt);CHKERRQ(ierr);
9811ae3d29cSBarry Smith     ierr = PetscMalloc(sizeof(char*),&amsopt->data);CHKERRQ(ierr);
982a297a907SKarl Rupp 
9831ae3d29cSBarry Smith     *(const char**)amsopt->data = defaultv;
984a297a907SKarl Rupp 
9851ae3d29cSBarry Smith     amsopt->list  = list;
9861ae3d29cSBarry Smith     amsopt->nlist = ntext;
9871ae3d29cSBarry Smith   }
98853acd3b1SBarry Smith   ierr = PetscOptionsGetEList(PetscOptionsObject.prefix,opt,list,ntext,value,set);CHKERRQ(ierr);
98961b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
99053acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%s> (choose one of)",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,defaultv);CHKERRQ(ierr);
99153acd3b1SBarry Smith     for (i=0; i<ntext; i++) {
99253acd3b1SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm," %s",list[i]);CHKERRQ(ierr);
99353acd3b1SBarry Smith     }
9942aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm," (%s)\n",ManSection(man));CHKERRQ(ierr);
99553acd3b1SBarry Smith   }
99653acd3b1SBarry Smith   PetscFunctionReturn(0);
99753acd3b1SBarry Smith }
99853acd3b1SBarry Smith 
99953acd3b1SBarry Smith #undef __FUNCT__
1000acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolGroupBegin"
100153acd3b1SBarry Smith /*@C
1002acfcf0e5SJed Brown      PetscOptionsBoolGroupBegin - First in a series of logical queries on the options database for
1003d5649816SBarry Smith        which at most a single value can be true.
100453acd3b1SBarry Smith 
10053f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
100653acd3b1SBarry Smith 
100753acd3b1SBarry Smith    Input Parameters:
100853acd3b1SBarry Smith +  opt - option name
100953acd3b1SBarry Smith .  text - short string that describes the option
101053acd3b1SBarry Smith -  man - manual page with additional information on option
101153acd3b1SBarry Smith 
101253acd3b1SBarry Smith    Output Parameter:
101353acd3b1SBarry Smith .  flg - whether that option was set or not
101453acd3b1SBarry Smith 
101553acd3b1SBarry Smith    Level: intermediate
101653acd3b1SBarry Smith 
101753acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
101853acd3b1SBarry Smith 
1019acfcf0e5SJed Brown    Must be followed by 0 or more PetscOptionsBoolGroup()s and PetscOptionsBoolGroupEnd()
102053acd3b1SBarry Smith 
102153acd3b1SBarry Smith     Concepts: options database^logical group
102253acd3b1SBarry Smith 
102353acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1024acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
102553acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
102653acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1027acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
102853acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
102953acd3b1SBarry Smith @*/
10307087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolGroupBegin(const char opt[],const char text[],const char man[],PetscBool  *flg)
103153acd3b1SBarry Smith {
103253acd3b1SBarry Smith   PetscErrorCode ierr;
10331ae3d29cSBarry Smith   PetscOptions   amsopt;
103453acd3b1SBarry Smith 
103553acd3b1SBarry Smith   PetscFunctionBegin;
10361ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
10371ae3d29cSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_LOGICAL,&amsopt);CHKERRQ(ierr);
1038ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1039a297a907SKarl Rupp 
1040ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
10411ae3d29cSBarry Smith   }
104268b16fdaSBarry Smith   *flg = PETSC_FALSE;
1043*0298fd71SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,NULL);CHKERRQ(ierr);
104461b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
104553acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  Pick at most one of -------------\n");CHKERRQ(ierr);
10462aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"    -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
104753acd3b1SBarry Smith   }
104853acd3b1SBarry Smith   PetscFunctionReturn(0);
104953acd3b1SBarry Smith }
105053acd3b1SBarry Smith 
105153acd3b1SBarry Smith #undef __FUNCT__
1052acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolGroup"
105353acd3b1SBarry Smith /*@C
1054acfcf0e5SJed Brown      PetscOptionsBoolGroup - One in a series of logical queries on the options database for
1055d5649816SBarry Smith        which at most a single value can be true.
105653acd3b1SBarry Smith 
10573f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
105853acd3b1SBarry Smith 
105953acd3b1SBarry Smith    Input Parameters:
106053acd3b1SBarry Smith +  opt - option name
106153acd3b1SBarry Smith .  text - short string that describes the option
106253acd3b1SBarry Smith -  man - manual page with additional information on option
106353acd3b1SBarry Smith 
106453acd3b1SBarry Smith    Output Parameter:
106553acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
106653acd3b1SBarry Smith 
106753acd3b1SBarry Smith    Level: intermediate
106853acd3b1SBarry Smith 
106953acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
107053acd3b1SBarry Smith 
1071acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin() and preceded a PetscOptionsBoolGroupEnd()
107253acd3b1SBarry Smith 
107353acd3b1SBarry Smith     Concepts: options database^logical group
107453acd3b1SBarry Smith 
107553acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1076acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
107753acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
107853acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1079acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
108053acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
108153acd3b1SBarry Smith @*/
10827087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolGroup(const char opt[],const char text[],const char man[],PetscBool  *flg)
108353acd3b1SBarry Smith {
108453acd3b1SBarry Smith   PetscErrorCode ierr;
10851ae3d29cSBarry Smith   PetscOptions   amsopt;
108653acd3b1SBarry Smith 
108753acd3b1SBarry Smith   PetscFunctionBegin;
10881ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
10891ae3d29cSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_LOGICAL,&amsopt);CHKERRQ(ierr);
1090ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1091a297a907SKarl Rupp 
1092ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
10931ae3d29cSBarry Smith   }
109417326d04SJed Brown   *flg = PETSC_FALSE;
1095*0298fd71SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,NULL);CHKERRQ(ierr);
109661b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
10972aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"    -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
109853acd3b1SBarry Smith   }
109953acd3b1SBarry Smith   PetscFunctionReturn(0);
110053acd3b1SBarry Smith }
110153acd3b1SBarry Smith 
110253acd3b1SBarry Smith #undef __FUNCT__
1103acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolGroupEnd"
110453acd3b1SBarry Smith /*@C
1105acfcf0e5SJed Brown      PetscOptionsBoolGroupEnd - Last in a series of logical queries on the options database for
1106d5649816SBarry Smith        which at most a single value can be true.
110753acd3b1SBarry Smith 
11083f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
110953acd3b1SBarry Smith 
111053acd3b1SBarry Smith    Input Parameters:
111153acd3b1SBarry Smith +  opt - option name
111253acd3b1SBarry Smith .  text - short string that describes the option
111353acd3b1SBarry Smith -  man - manual page with additional information on option
111453acd3b1SBarry Smith 
111553acd3b1SBarry Smith    Output Parameter:
111653acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
111753acd3b1SBarry Smith 
111853acd3b1SBarry Smith    Level: intermediate
111953acd3b1SBarry Smith 
112053acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
112153acd3b1SBarry Smith 
1122acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin()
112353acd3b1SBarry Smith 
112453acd3b1SBarry Smith     Concepts: options database^logical group
112553acd3b1SBarry Smith 
112653acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1127acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
112853acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
112953acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1130acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
113153acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
113253acd3b1SBarry Smith @*/
11337087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolGroupEnd(const char opt[],const char text[],const char man[],PetscBool  *flg)
113453acd3b1SBarry Smith {
113553acd3b1SBarry Smith   PetscErrorCode ierr;
11361ae3d29cSBarry Smith   PetscOptions   amsopt;
113753acd3b1SBarry Smith 
113853acd3b1SBarry Smith   PetscFunctionBegin;
11391ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
11401ae3d29cSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_LOGICAL,&amsopt);CHKERRQ(ierr);
1141ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1142a297a907SKarl Rupp 
1143ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
11441ae3d29cSBarry Smith   }
114517326d04SJed Brown   *flg = PETSC_FALSE;
1146*0298fd71SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,NULL);CHKERRQ(ierr);
114761b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
11482aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"    -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
114953acd3b1SBarry Smith   }
115053acd3b1SBarry Smith   PetscFunctionReturn(0);
115153acd3b1SBarry Smith }
115253acd3b1SBarry Smith 
115353acd3b1SBarry Smith #undef __FUNCT__
1154acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBool"
115553acd3b1SBarry Smith /*@C
1156acfcf0e5SJed Brown    PetscOptionsBool - Determines if a particular option is in the database with a true or false
115753acd3b1SBarry Smith 
11583f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
115953acd3b1SBarry Smith 
116053acd3b1SBarry Smith    Input Parameters:
116153acd3b1SBarry Smith +  opt - option name
116253acd3b1SBarry Smith .  text - short string that describes the option
116353acd3b1SBarry Smith -  man - manual page with additional information on option
116453acd3b1SBarry Smith 
116553acd3b1SBarry Smith    Output Parameter:
116653acd3b1SBarry Smith .  flg - PETSC_TRUE or PETSC_FALSE
116753acd3b1SBarry Smith .  set - PETSC_TRUE if found, else PETSC_FALSE
116853acd3b1SBarry Smith 
116953acd3b1SBarry Smith    Level: beginner
117053acd3b1SBarry Smith 
117153acd3b1SBarry Smith    Concepts: options database^logical
117253acd3b1SBarry Smith 
117353acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
117453acd3b1SBarry Smith 
117553acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
1176acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
1177acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
117853acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
117953acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1180acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
118153acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
118253acd3b1SBarry Smith @*/
11837087cfbeSBarry Smith PetscErrorCode  PetscOptionsBool(const char opt[],const char text[],const char man[],PetscBool deflt,PetscBool  *flg,PetscBool  *set)
118453acd3b1SBarry Smith {
118553acd3b1SBarry Smith   PetscErrorCode ierr;
1186ace3abfcSBarry Smith   PetscBool      iset;
1187aee2cecaSBarry Smith   PetscOptions   amsopt;
118853acd3b1SBarry Smith 
118953acd3b1SBarry Smith   PetscFunctionBegin;
1190b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
1191aee2cecaSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_LOGICAL,&amsopt);CHKERRQ(ierr);
1192ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1193a297a907SKarl Rupp 
1194ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = deflt;
1195af6d86caSBarry Smith   }
1196acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,&iset);CHKERRQ(ierr);
119753acd3b1SBarry Smith   if (!iset) {
119853acd3b1SBarry Smith     if (flg) *flg = deflt;
119953acd3b1SBarry Smith   }
120053acd3b1SBarry Smith   if (set) *set = iset;
120161b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
1202ace3abfcSBarry Smith     const char *v = PetscBools[deflt];
12032aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s: <%s> %s (%s)\n",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,v,text,ManSection(man));CHKERRQ(ierr);
120453acd3b1SBarry Smith   }
120553acd3b1SBarry Smith   PetscFunctionReturn(0);
120653acd3b1SBarry Smith }
120753acd3b1SBarry Smith 
120853acd3b1SBarry Smith #undef __FUNCT__
120953acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsRealArray"
121053acd3b1SBarry Smith /*@C
121153acd3b1SBarry Smith    PetscOptionsRealArray - Gets an array of double values for a particular
121253acd3b1SBarry Smith    option in the database. The values must be separated with commas with
121353acd3b1SBarry Smith    no intervening spaces.
121453acd3b1SBarry Smith 
12153f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
121653acd3b1SBarry Smith 
121753acd3b1SBarry Smith    Input Parameters:
121853acd3b1SBarry Smith +  opt - the option one is seeking
121953acd3b1SBarry Smith .  text - short string describing option
122053acd3b1SBarry Smith .  man - manual page for option
122153acd3b1SBarry Smith -  nmax - maximum number of values
122253acd3b1SBarry Smith 
122353acd3b1SBarry Smith    Output Parameter:
122453acd3b1SBarry Smith +  value - location to copy values
122553acd3b1SBarry Smith .  nmax - actual number of values found
122653acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
122753acd3b1SBarry Smith 
122853acd3b1SBarry Smith    Level: beginner
122953acd3b1SBarry Smith 
123053acd3b1SBarry Smith    Notes:
123153acd3b1SBarry Smith    The user should pass in an array of doubles
123253acd3b1SBarry Smith 
123353acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
123453acd3b1SBarry Smith 
123553acd3b1SBarry Smith    Concepts: options database^array of strings
123653acd3b1SBarry Smith 
123753acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1238acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
123953acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
124053acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1241acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
124253acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
124353acd3b1SBarry Smith @*/
12447087cfbeSBarry Smith PetscErrorCode  PetscOptionsRealArray(const char opt[],const char text[],const char man[],PetscReal value[],PetscInt *n,PetscBool  *set)
124553acd3b1SBarry Smith {
124653acd3b1SBarry Smith   PetscErrorCode ierr;
124753acd3b1SBarry Smith   PetscInt       i;
1248e26ddf31SBarry Smith   PetscOptions   amsopt;
124953acd3b1SBarry Smith 
125053acd3b1SBarry Smith   PetscFunctionBegin;
1251b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
1252e26ddf31SBarry Smith     PetscReal *vals;
1253e26ddf31SBarry Smith 
1254e26ddf31SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_REAL_ARRAY,&amsopt);CHKERRQ(ierr);
1255e26ddf31SBarry Smith     ierr = PetscMalloc((*n)*sizeof(PetscReal),&amsopt->data);CHKERRQ(ierr);
1256e26ddf31SBarry Smith     vals = (PetscReal*)amsopt->data;
1257e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1258e26ddf31SBarry Smith     amsopt->arraylength = *n;
1259e26ddf31SBarry Smith   }
126053acd3b1SBarry Smith   ierr = PetscOptionsGetRealArray(PetscOptionsObject.prefix,opt,value,n,set);CHKERRQ(ierr);
126161b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
1262a83599f4SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%G",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,value[0]);CHKERRQ(ierr);
126353acd3b1SBarry Smith     for (i=1; i<*n; i++) {
1264a83599f4SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,",%G",value[i]);CHKERRQ(ierr);
126553acd3b1SBarry Smith     }
12662aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
126753acd3b1SBarry Smith   }
126853acd3b1SBarry Smith   PetscFunctionReturn(0);
126953acd3b1SBarry Smith }
127053acd3b1SBarry Smith 
127153acd3b1SBarry Smith 
127253acd3b1SBarry Smith #undef __FUNCT__
127353acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsIntArray"
127453acd3b1SBarry Smith /*@C
127553acd3b1SBarry Smith    PetscOptionsIntArray - Gets an array of integers for a particular
1276b32a342fSShri Abhyankar    option in the database.
127753acd3b1SBarry Smith 
12783f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
127953acd3b1SBarry Smith 
128053acd3b1SBarry Smith    Input Parameters:
128153acd3b1SBarry Smith +  opt - the option one is seeking
128253acd3b1SBarry Smith .  text - short string describing option
128353acd3b1SBarry Smith .  man - manual page for option
1284f8a50e2bSBarry Smith -  n - maximum number of values
128553acd3b1SBarry Smith 
128653acd3b1SBarry Smith    Output Parameter:
128753acd3b1SBarry Smith +  value - location to copy values
1288f8a50e2bSBarry Smith .  n - actual number of values found
128953acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
129053acd3b1SBarry Smith 
129153acd3b1SBarry Smith    Level: beginner
129253acd3b1SBarry Smith 
129353acd3b1SBarry Smith    Notes:
1294b32a342fSShri Abhyankar    The array can be passed as
1295b32a342fSShri Abhyankar    a comma seperated list:                                 0,1,2,3,4,5,6,7
12960fd488f5SShri Abhyankar    a range (start-end+1):                                  0-8
12970fd488f5SShri Abhyankar    a range with given increment (start-end+1:inc):         0-7:2
12980fd488f5SShri Abhyankar    a combination of values and ranges seperated by commas: 0,1-8,8-15:2
1299b32a342fSShri Abhyankar 
1300b32a342fSShri Abhyankar    There must be no intervening spaces between the values.
130153acd3b1SBarry Smith 
130253acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
130353acd3b1SBarry Smith 
1304b32a342fSShri Abhyankar    Concepts: options database^array of ints
130553acd3b1SBarry Smith 
130653acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1307acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
130853acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
130953acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1310acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
131153acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList(), PetscOptionsRealArray()
131253acd3b1SBarry Smith @*/
13137087cfbeSBarry Smith PetscErrorCode  PetscOptionsIntArray(const char opt[],const char text[],const char man[],PetscInt value[],PetscInt *n,PetscBool  *set)
131453acd3b1SBarry Smith {
131553acd3b1SBarry Smith   PetscErrorCode ierr;
131653acd3b1SBarry Smith   PetscInt       i;
1317e26ddf31SBarry Smith   PetscOptions   amsopt;
131853acd3b1SBarry Smith 
131953acd3b1SBarry Smith   PetscFunctionBegin;
1320b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
1321e26ddf31SBarry Smith     PetscInt *vals;
1322e26ddf31SBarry Smith 
1323e26ddf31SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_INT_ARRAY,&amsopt);CHKERRQ(ierr);
1324e26ddf31SBarry Smith     ierr = PetscMalloc((*n)*sizeof(PetscInt),&amsopt->data);CHKERRQ(ierr);
1325e26ddf31SBarry Smith     vals = (PetscInt*)amsopt->data;
1326e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1327e26ddf31SBarry Smith     amsopt->arraylength = *n;
1328e26ddf31SBarry Smith   }
132953acd3b1SBarry Smith   ierr = PetscOptionsGetIntArray(PetscOptionsObject.prefix,opt,value,n,set);CHKERRQ(ierr);
133061b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
133153acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%d",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,value[0]);CHKERRQ(ierr);
133253acd3b1SBarry Smith     for (i=1; i<*n; i++) {
133353acd3b1SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,",%d",value[i]);CHKERRQ(ierr);
133453acd3b1SBarry Smith     }
13352aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
133653acd3b1SBarry Smith   }
133753acd3b1SBarry Smith   PetscFunctionReturn(0);
133853acd3b1SBarry Smith }
133953acd3b1SBarry Smith 
134053acd3b1SBarry Smith #undef __FUNCT__
134153acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsStringArray"
134253acd3b1SBarry Smith /*@C
134353acd3b1SBarry Smith    PetscOptionsStringArray - Gets an array of string values for a particular
134453acd3b1SBarry Smith    option in the database. The values must be separated with commas with
134553acd3b1SBarry Smith    no intervening spaces.
134653acd3b1SBarry Smith 
13473f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
134853acd3b1SBarry Smith 
134953acd3b1SBarry Smith    Input Parameters:
135053acd3b1SBarry Smith +  opt - the option one is seeking
135153acd3b1SBarry Smith .  text - short string describing option
135253acd3b1SBarry Smith .  man - manual page for option
135353acd3b1SBarry Smith -  nmax - maximum number of strings
135453acd3b1SBarry Smith 
135553acd3b1SBarry Smith    Output Parameter:
135653acd3b1SBarry Smith +  value - location to copy strings
135753acd3b1SBarry Smith .  nmax - actual number of strings found
135853acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
135953acd3b1SBarry Smith 
136053acd3b1SBarry Smith    Level: beginner
136153acd3b1SBarry Smith 
136253acd3b1SBarry Smith    Notes:
136353acd3b1SBarry Smith    The user should pass in an array of pointers to char, to hold all the
136453acd3b1SBarry Smith    strings returned by this function.
136553acd3b1SBarry Smith 
136653acd3b1SBarry Smith    The user is responsible for deallocating the strings that are
136753acd3b1SBarry Smith    returned. The Fortran interface for this routine is not supported.
136853acd3b1SBarry Smith 
136953acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
137053acd3b1SBarry Smith 
137153acd3b1SBarry Smith    Concepts: options database^array of strings
137253acd3b1SBarry Smith 
137353acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1374acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
137553acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
137653acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1377acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
137853acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
137953acd3b1SBarry Smith @*/
13807087cfbeSBarry Smith PetscErrorCode  PetscOptionsStringArray(const char opt[],const char text[],const char man[],char *value[],PetscInt *nmax,PetscBool  *set)
138153acd3b1SBarry Smith {
138253acd3b1SBarry Smith   PetscErrorCode ierr;
13831ae3d29cSBarry Smith   PetscOptions   amsopt;
138453acd3b1SBarry Smith 
138553acd3b1SBarry Smith   PetscFunctionBegin;
13861ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
13871ae3d29cSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_STRING_ARRAY,&amsopt);CHKERRQ(ierr);
13881ae3d29cSBarry Smith     ierr = PetscMalloc((*nmax)*sizeof(char*),&amsopt->data);CHKERRQ(ierr);
1389a297a907SKarl Rupp 
13901ae3d29cSBarry Smith     amsopt->arraylength = *nmax;
13911ae3d29cSBarry Smith   }
139253acd3b1SBarry Smith   ierr = PetscOptionsGetStringArray(PetscOptionsObject.prefix,opt,value,nmax,set);CHKERRQ(ierr);
139361b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
13942aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <string1,string2,...>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
139553acd3b1SBarry Smith   }
139653acd3b1SBarry Smith   PetscFunctionReturn(0);
139753acd3b1SBarry Smith }
139853acd3b1SBarry Smith 
1399e2446a98SMatthew Knepley #undef __FUNCT__
1400acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolArray"
1401e2446a98SMatthew Knepley /*@C
1402acfcf0e5SJed Brown    PetscOptionsBoolArray - Gets an array of logical values (true or false) for a particular
1403e2446a98SMatthew Knepley    option in the database. The values must be separated with commas with
1404e2446a98SMatthew Knepley    no intervening spaces.
1405e2446a98SMatthew Knepley 
14063f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
1407e2446a98SMatthew Knepley 
1408e2446a98SMatthew Knepley    Input Parameters:
1409e2446a98SMatthew Knepley +  opt - the option one is seeking
1410e2446a98SMatthew Knepley .  text - short string describing option
1411e2446a98SMatthew Knepley .  man - manual page for option
1412e2446a98SMatthew Knepley -  nmax - maximum number of values
1413e2446a98SMatthew Knepley 
1414e2446a98SMatthew Knepley    Output Parameter:
1415e2446a98SMatthew Knepley +  value - location to copy values
1416e2446a98SMatthew Knepley .  nmax - actual number of values found
1417e2446a98SMatthew Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
1418e2446a98SMatthew Knepley 
1419e2446a98SMatthew Knepley    Level: beginner
1420e2446a98SMatthew Knepley 
1421e2446a98SMatthew Knepley    Notes:
1422e2446a98SMatthew Knepley    The user should pass in an array of doubles
1423e2446a98SMatthew Knepley 
1424e2446a98SMatthew Knepley    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
1425e2446a98SMatthew Knepley 
1426e2446a98SMatthew Knepley    Concepts: options database^array of strings
1427e2446a98SMatthew Knepley 
1428e2446a98SMatthew Knepley .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1429acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
1430e2446a98SMatthew Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1431e2446a98SMatthew Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1432acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1433e2446a98SMatthew Knepley           PetscOptionsList(), PetscOptionsEList()
1434e2446a98SMatthew Knepley @*/
14357087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolArray(const char opt[],const char text[],const char man[],PetscBool value[],PetscInt *n,PetscBool *set)
1436e2446a98SMatthew Knepley {
1437e2446a98SMatthew Knepley   PetscErrorCode ierr;
1438e2446a98SMatthew Knepley   PetscInt       i;
14391ae3d29cSBarry Smith   PetscOptions   amsopt;
1440e2446a98SMatthew Knepley 
1441e2446a98SMatthew Knepley   PetscFunctionBegin;
14421ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
1443ace3abfcSBarry Smith     PetscBool *vals;
14441ae3d29cSBarry Smith 
14451ae3d29cSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_LOGICAL_ARRAY,&amsopt);CHKERRQ(ierr);
1446ace3abfcSBarry Smith     ierr = PetscMalloc((*n)*sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1447ace3abfcSBarry Smith     vals = (PetscBool*)amsopt->data;
14481ae3d29cSBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
14491ae3d29cSBarry Smith     amsopt->arraylength = *n;
14501ae3d29cSBarry Smith   }
1451acfcf0e5SJed Brown   ierr = PetscOptionsGetBoolArray(PetscOptionsObject.prefix,opt,value,n,set);CHKERRQ(ierr);
1452e2446a98SMatthew Knepley   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
1453e2446a98SMatthew Knepley     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%d",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,value[0]);CHKERRQ(ierr);
1454e2446a98SMatthew Knepley     for (i=1; i<*n; i++) {
1455e2446a98SMatthew Knepley       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,",%d",value[i]);CHKERRQ(ierr);
1456e2446a98SMatthew Knepley     }
14572aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
1458e2446a98SMatthew Knepley   }
1459e2446a98SMatthew Knepley   PetscFunctionReturn(0);
1460e2446a98SMatthew Knepley }
1461e2446a98SMatthew Knepley 
14628cc676e6SMatthew G Knepley #undef __FUNCT__
14638cc676e6SMatthew G Knepley #define __FUNCT__ "PetscOptionsViewer"
14648cc676e6SMatthew G Knepley /*@C
14658cc676e6SMatthew G Knepley    PetscOptionsInt - Gets a viewer appropriate for the type indicated by the user
14668cc676e6SMatthew G Knepley 
14678cc676e6SMatthew G Knepley    Logically Collective on the communicator passed in PetscOptionsBegin()
14688cc676e6SMatthew G Knepley 
14698cc676e6SMatthew G Knepley    Input Parameters:
14708cc676e6SMatthew G Knepley +  opt - option name
14718cc676e6SMatthew G Knepley .  text - short string that describes the option
14728cc676e6SMatthew G Knepley -  man - manual page with additional information on option
14738cc676e6SMatthew G Knepley 
14748cc676e6SMatthew G Knepley    Output Parameter:
14758cc676e6SMatthew G Knepley +  viewer - the viewer
14768cc676e6SMatthew G Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
14778cc676e6SMatthew G Knepley 
14788cc676e6SMatthew G Knepley    Level: beginner
14798cc676e6SMatthew G Knepley 
14808cc676e6SMatthew G Knepley    Concepts: options database^has int
14818cc676e6SMatthew G Knepley 
14828cc676e6SMatthew G Knepley    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
14838cc676e6SMatthew G Knepley If no value is provided ascii:stdout is used
14848cc676e6SMatthew 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
14858cc676e6SMatthew G Knepley $                                     about the object to standard out
14868cc676e6SMatthew G Knepley $       binary[:filename]   defaults to binaryoutput
14878cc676e6SMatthew G Knepley $       draw
14888cc676e6SMatthew G Knepley $       socket[:port]    defaults to the standard output port
14898cc676e6SMatthew G Knepley 
1490cffb1e40SBarry Smith    Use PetscRestoreViewerDestroy() after using the viewer, otherwise a memory leak will occur
14918cc676e6SMatthew G Knepley 
14928cc676e6SMatthew G Knepley .seealso: PetscOptionsGetViewer(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
14938cc676e6SMatthew G Knepley           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
14948cc676e6SMatthew G Knepley           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
14958cc676e6SMatthew G Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
14968cc676e6SMatthew G Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
14978cc676e6SMatthew G Knepley           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
14988cc676e6SMatthew G Knepley           PetscOptionsList(), PetscOptionsEList()
14998cc676e6SMatthew G Knepley @*/
1500cffb1e40SBarry Smith PetscErrorCode  PetscOptionsViewer(const char opt[],const char text[],const char man[],PetscViewer *viewer,PetscViewerFormat *format,PetscBool  *set)
15018cc676e6SMatthew G Knepley {
15028cc676e6SMatthew G Knepley   PetscErrorCode ierr;
15038cc676e6SMatthew G Knepley   PetscOptions   amsopt;
15048cc676e6SMatthew G Knepley 
15058cc676e6SMatthew G Knepley   PetscFunctionBegin;
15068cc676e6SMatthew G Knepley   if (!PetscOptionsPublishCount) {
15078cc676e6SMatthew G Knepley     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
15088cc676e6SMatthew G Knepley     ierr = PetscMalloc(sizeof(PetscInt),&amsopt->data);CHKERRQ(ierr);
1509a297a907SKarl Rupp 
15108cc676e6SMatthew G Knepley     *(const char**)amsopt->data = "";
15118cc676e6SMatthew G Knepley   }
1512cffb1e40SBarry Smith   ierr = PetscOptionsGetViewer(PetscOptionsObject.comm,PetscOptionsObject.prefix,opt,viewer,format,set);CHKERRQ(ierr);
15138cc676e6SMatthew G Knepley   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
15148cc676e6SMatthew G Knepley     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%s>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,"",text,ManSection(man));CHKERRQ(ierr);
15158cc676e6SMatthew G Knepley   }
15168cc676e6SMatthew G Knepley   PetscFunctionReturn(0);
15178cc676e6SMatthew G Knepley }
15188cc676e6SMatthew G Knepley 
151953acd3b1SBarry Smith 
152053acd3b1SBarry Smith #undef __FUNCT__
152153acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsHead"
152253acd3b1SBarry Smith /*@C
1523b52f573bSBarry Smith      PetscOptionsHead - Puts a heading before listing any more published options. Used, for example,
152453acd3b1SBarry Smith             in KSPSetFromOptions_GMRES().
152553acd3b1SBarry Smith 
15263f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
152753acd3b1SBarry Smith 
152853acd3b1SBarry Smith    Input Parameter:
152953acd3b1SBarry Smith .   head - the heading text
153053acd3b1SBarry Smith 
153153acd3b1SBarry Smith 
153253acd3b1SBarry Smith    Level: intermediate
153353acd3b1SBarry Smith 
153453acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
153553acd3b1SBarry Smith 
1536b52f573bSBarry Smith           Can be followed by a call to PetscOptionsTail() in the same function.
153753acd3b1SBarry Smith 
153853acd3b1SBarry Smith    Concepts: options database^subheading
153953acd3b1SBarry Smith 
154053acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1541acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
154253acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
154353acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1544acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
154553acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
154653acd3b1SBarry Smith @*/
15477087cfbeSBarry Smith PetscErrorCode  PetscOptionsHead(const char head[])
154853acd3b1SBarry Smith {
154953acd3b1SBarry Smith   PetscErrorCode ierr;
155053acd3b1SBarry Smith 
155153acd3b1SBarry Smith   PetscFunctionBegin;
155261b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
155353acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  %s\n",head);CHKERRQ(ierr);
155453acd3b1SBarry Smith   }
155553acd3b1SBarry Smith   PetscFunctionReturn(0);
155653acd3b1SBarry Smith }
155753acd3b1SBarry Smith 
155853acd3b1SBarry Smith 
155953acd3b1SBarry Smith 
156053acd3b1SBarry Smith 
156153acd3b1SBarry Smith 
156253acd3b1SBarry Smith 
1563