xref: /petsc/src/sys/objects/aoptions.c (revision 9de0f6ec9b026c9f9bc5a81a4cedd1340deba9ac)
17d0a6c19SBarry Smith 
21a1499c8SBarry Smith 
353acd3b1SBarry Smith /*
43fc1eb6aSBarry Smith    Implements the higher-level options database querying methods. These are self-documenting and can attach at runtime to
53fc1eb6aSBarry Smith    GUI code to display the options and get values from the users.
653acd3b1SBarry Smith 
753acd3b1SBarry Smith */
853acd3b1SBarry Smith 
9af0996ceSBarry Smith #include <petsc/private/petscimpl.h>        /*I  "petscsys.h"   I*/
10665c2dedSJed Brown #include <petscviewer.h>
1153acd3b1SBarry Smith 
122aa6d131SJed Brown #define ManSection(str) ((str) ? (str) : "None")
132aa6d131SJed Brown 
1453acd3b1SBarry Smith /*
1553acd3b1SBarry Smith     Keep a linked list of options that have been posted and we are waiting for
163fc1eb6aSBarry Smith    user selection. See the manual page for PetscOptionsBegin()
1753acd3b1SBarry Smith 
1853acd3b1SBarry Smith     Eventually we'll attach this beast to a MPI_Comm
1953acd3b1SBarry Smith */
20e55864a3SBarry Smith 
2153acd3b1SBarry Smith 
2253acd3b1SBarry Smith #undef __FUNCT__
2353acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsBegin_Private"
2453acd3b1SBarry Smith /*
2553acd3b1SBarry Smith     Handles setting up the data structure in a call to PetscOptionsBegin()
2653acd3b1SBarry Smith */
274416b707SBarry Smith PetscErrorCode PetscOptionsBegin_Private(PetscOptionItems *PetscOptionsObject,MPI_Comm comm,const char prefix[],const char title[],const char mansec[])
2853acd3b1SBarry Smith {
2953acd3b1SBarry Smith   PetscErrorCode ierr;
3053acd3b1SBarry Smith 
3153acd3b1SBarry Smith   PetscFunctionBegin;
32*9de0f6ecSBarry Smith   if (!PetscOptionsHelpPrintedSingleton) {
33*9de0f6ecSBarry Smith     ierr = PetscOptionsHelpPrintedCreate(&PetscOptionsHelpPrintedSingleton);CHKERRQ(ierr);
34*9de0f6ecSBarry Smith   }
35*9de0f6ecSBarry Smith   ierr = PetscOptionsHelpPrintedCheck(PetscOptionsHelpPrintedSingleton,prefix,title,&PetscOptionsObject->alreadyprinted);CHKERRQ(ierr);
36e55864a3SBarry Smith   PetscOptionsObject->next          = 0;
37e55864a3SBarry Smith   PetscOptionsObject->comm          = comm;
38e55864a3SBarry Smith   PetscOptionsObject->changedmethod = PETSC_FALSE;
39a297a907SKarl Rupp 
40e55864a3SBarry Smith   ierr = PetscStrallocpy(prefix,&PetscOptionsObject->prefix);CHKERRQ(ierr);
41e55864a3SBarry Smith   ierr = PetscStrallocpy(title,&PetscOptionsObject->title);CHKERRQ(ierr);
4253acd3b1SBarry Smith 
43c5929fdfSBarry Smith   ierr = PetscOptionsHasName(PetscOptionsObject->options,NULL,"-help",&PetscOptionsObject->printhelp);CHKERRQ(ierr);
44e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1) {
45e55864a3SBarry Smith     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 */
574416b707SBarry Smith PetscErrorCode PetscObjectOptionsBegin_Private(PetscOptionItems *PetscOptionsObject,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);
65e55864a3SBarry Smith   PetscOptionsObject->object         = obj;
66e55864a3SBarry Smith   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   }
74e55864a3SBarry Smith   ierr = PetscOptionsBegin_Private(PetscOptionsObject,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__
824416b707SBarry Smith #define __FUNCT__ "PetscOptionCreate_Private"
834416b707SBarry Smith static int PetscOptionItemCreate_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscOptionType t,PetscOptionItem *amsopt)
8453acd3b1SBarry Smith {
8553acd3b1SBarry Smith   int             ierr;
864416b707SBarry Smith   PetscOptionItem 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 
93b00a9115SJed Brown   ierr            = PetscNew(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 
103e55864a3SBarry Smith   if (!PetscOptionsObject->next) PetscOptionsObject->next = *amsopt;
104a297a907SKarl Rupp   else {
105e55864a3SBarry 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 
151ead66b60SBarry Smith #undef __FUNCT__
152ead66b60SBarry Smith #define __FUNCT__ "PetscStrdup"
1535b02f95dSBarry Smith /*
1545b02f95dSBarry Smith     This is needed because certain strings may be freed by SAWs, hence we cannot use PetscStrallocpy()
1555b02f95dSBarry Smith */
1565b02f95dSBarry Smith static PetscErrorCode  PetscStrdup(const char s[],char *t[])
1575b02f95dSBarry Smith {
1585b02f95dSBarry Smith   PetscErrorCode ierr;
1595b02f95dSBarry Smith   size_t         len;
1605b02f95dSBarry Smith   char           *tmp = 0;
1615b02f95dSBarry Smith 
1625b02f95dSBarry Smith   PetscFunctionBegin;
1635b02f95dSBarry Smith   if (s) {
1645b02f95dSBarry Smith     ierr = PetscStrlen(s,&len);CHKERRQ(ierr);
165f416af30SBarry Smith     tmp = (char*) malloc((len+1)*sizeof(char));
1665b02f95dSBarry Smith     if (!tmp) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"No memory to duplicate string");
1675b02f95dSBarry Smith     ierr = PetscStrcpy(tmp,s);CHKERRQ(ierr);
1685b02f95dSBarry Smith   }
1695b02f95dSBarry Smith   *t = tmp;
1705b02f95dSBarry Smith   PetscFunctionReturn(0);
1715b02f95dSBarry Smith }
1725b02f95dSBarry Smith 
1735b02f95dSBarry Smith 
174aee2cecaSBarry Smith #undef __FUNCT__
175aee2cecaSBarry Smith #define __FUNCT__ "PetscOptionsGetFromTextInput"
176aee2cecaSBarry Smith /*
1773cc1e11dSBarry Smith     PetscOptionsGetFromTextInput - Presents all the PETSc Options processed by the program so the user may change them at runtime
178aee2cecaSBarry Smith 
179aee2cecaSBarry Smith     Notes: this isn't really practical, it is just to demonstrate the principle
180aee2cecaSBarry Smith 
1817781c08eSBarry Smith     A carriage return indicates no change from the default; but this like -ksp_monitor <stdout>  the default is actually not stdout the default
1827781c08eSBarry Smith     is to do nothing so to get it to use stdout you need to type stdout. This is kind of bug?
1837781c08eSBarry Smith 
184aee2cecaSBarry Smith     Bugs:
1857781c08eSBarry Smith +    All processes must traverse through the exact same set of option queries due to the call to PetscScanString()
1863cc1e11dSBarry Smith .    Internal strings have arbitrary length and string copies are not checked that they fit into string space
187aee2cecaSBarry Smith -    Only works for PetscInt == int, PetscReal == double etc
188aee2cecaSBarry Smith 
1893cc1e11dSBarry Smith     Developer Notes: Normally the GUI that presents the options the user and retrieves the values would be running in a different
1903cc1e11dSBarry Smith      address space and communicating with the PETSc program
1913cc1e11dSBarry Smith 
192aee2cecaSBarry Smith */
1934416b707SBarry Smith PetscErrorCode PetscOptionsGetFromTextInput(PetscOptionItems *PetscOptionsObject)
1946356e834SBarry Smith {
1956356e834SBarry Smith   PetscErrorCode  ierr;
1964416b707SBarry Smith   PetscOptionItem next = PetscOptionsObject->next;
1976356e834SBarry Smith   char            str[512];
1987781c08eSBarry Smith   PetscBool       bid;
199a4404d99SBarry Smith   PetscReal       ir,*valr;
200330cf3c9SBarry Smith   PetscInt        *vald;
201330cf3c9SBarry Smith   size_t          i;
2026356e834SBarry Smith 
2032a409bb0SBarry Smith   PetscFunctionBegin;
204e55864a3SBarry Smith   ierr = (*PetscPrintf)(PETSC_COMM_WORLD,"%s -------------------------------------------------\n",PetscOptionsObject->title);CHKERRQ(ierr);
2056356e834SBarry Smith   while (next) {
2066356e834SBarry Smith     switch (next->type) {
2076356e834SBarry Smith     case OPTION_HEAD:
2086356e834SBarry Smith       break;
209e26ddf31SBarry Smith     case OPTION_INT_ARRAY:
210e55864a3SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",next->option+1);CHKERRQ(ierr);
211e26ddf31SBarry Smith       vald = (PetscInt*) next->data;
212e26ddf31SBarry Smith       for (i=0; i<next->arraylength; i++) {
213e26ddf31SBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"%d",vald[i]);CHKERRQ(ierr);
214e26ddf31SBarry Smith         if (i < next->arraylength-1) {
215e26ddf31SBarry Smith           ierr = PetscPrintf(PETSC_COMM_WORLD,",");CHKERRQ(ierr);
216e26ddf31SBarry Smith         }
217e26ddf31SBarry Smith       }
218e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,">: %s (%s) ",next->text,next->man);CHKERRQ(ierr);
219e26ddf31SBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
220e26ddf31SBarry Smith       if (str[0]) {
221e26ddf31SBarry Smith         PetscToken token;
222e26ddf31SBarry Smith         PetscInt   n=0,nmax = next->arraylength,*dvalue = (PetscInt*)next->data,start,end;
223e26ddf31SBarry Smith         size_t     len;
224e26ddf31SBarry Smith         char       *value;
225ace3abfcSBarry Smith         PetscBool  foundrange;
226e26ddf31SBarry Smith 
227e26ddf31SBarry Smith         next->set = PETSC_TRUE;
228e26ddf31SBarry Smith         value     = str;
229e26ddf31SBarry Smith         ierr      = PetscTokenCreate(value,',',&token);CHKERRQ(ierr);
230e26ddf31SBarry Smith         ierr      = PetscTokenFind(token,&value);CHKERRQ(ierr);
231e26ddf31SBarry Smith         while (n < nmax) {
232e26ddf31SBarry Smith           if (!value) break;
233e26ddf31SBarry Smith 
234e26ddf31SBarry Smith           /* look for form  d-D where d and D are integers */
235e26ddf31SBarry Smith           foundrange = PETSC_FALSE;
236e26ddf31SBarry Smith           ierr       = PetscStrlen(value,&len);CHKERRQ(ierr);
237e26ddf31SBarry Smith           if (value[0] == '-') i=2;
238e26ddf31SBarry Smith           else i=1;
239330cf3c9SBarry Smith           for (;i<len; i++) {
240e26ddf31SBarry Smith             if (value[i] == '-') {
241e32f2f54SBarry Smith               if (i == len-1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Error in %D-th array entry %s\n",n,value);
242e26ddf31SBarry Smith               value[i] = 0;
243cfbddea1SSatish Balay               ierr     = PetscOptionsStringToInt(value,&start);CHKERRQ(ierr);
244cfbddea1SSatish Balay               ierr     = PetscOptionsStringToInt(value+i+1,&end);CHKERRQ(ierr);
245e32f2f54SBarry 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);
246e32f2f54SBarry 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);
247e26ddf31SBarry Smith               for (; start<end; start++) {
248e26ddf31SBarry Smith                 *dvalue = start; dvalue++;n++;
249e26ddf31SBarry Smith               }
250e26ddf31SBarry Smith               foundrange = PETSC_TRUE;
251e26ddf31SBarry Smith               break;
252e26ddf31SBarry Smith             }
253e26ddf31SBarry Smith           }
254e26ddf31SBarry Smith           if (!foundrange) {
255cfbddea1SSatish Balay             ierr = PetscOptionsStringToInt(value,dvalue);CHKERRQ(ierr);
256e26ddf31SBarry Smith             dvalue++;
257e26ddf31SBarry Smith             n++;
258e26ddf31SBarry Smith           }
259e26ddf31SBarry Smith           ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
260e26ddf31SBarry Smith         }
2618c74ee41SBarry Smith         ierr = PetscTokenDestroy(&token);CHKERRQ(ierr);
262e26ddf31SBarry Smith       }
263e26ddf31SBarry Smith       break;
264e26ddf31SBarry Smith     case OPTION_REAL_ARRAY:
265e55864a3SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",next->option+1);CHKERRQ(ierr);
266e26ddf31SBarry Smith       valr = (PetscReal*) next->data;
267e26ddf31SBarry Smith       for (i=0; i<next->arraylength; i++) {
268e26ddf31SBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"%g",valr[i]);CHKERRQ(ierr);
269e26ddf31SBarry Smith         if (i < next->arraylength-1) {
270e26ddf31SBarry Smith           ierr = PetscPrintf(PETSC_COMM_WORLD,",");CHKERRQ(ierr);
271e26ddf31SBarry Smith         }
272e26ddf31SBarry Smith       }
273e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,">: %s (%s) ",next->text,next->man);CHKERRQ(ierr);
274e26ddf31SBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
275e26ddf31SBarry Smith       if (str[0]) {
276e26ddf31SBarry Smith         PetscToken token;
277e26ddf31SBarry Smith         PetscInt   n = 0,nmax = next->arraylength;
278e26ddf31SBarry Smith         PetscReal  *dvalue = (PetscReal*)next->data;
279e26ddf31SBarry Smith         char       *value;
280e26ddf31SBarry Smith 
281e26ddf31SBarry Smith         next->set = PETSC_TRUE;
282e26ddf31SBarry Smith         value     = str;
283e26ddf31SBarry Smith         ierr      = PetscTokenCreate(value,',',&token);CHKERRQ(ierr);
284e26ddf31SBarry Smith         ierr      = PetscTokenFind(token,&value);CHKERRQ(ierr);
285e26ddf31SBarry Smith         while (n < nmax) {
286e26ddf31SBarry Smith           if (!value) break;
287cfbddea1SSatish Balay           ierr = PetscOptionsStringToReal(value,dvalue);CHKERRQ(ierr);
288e26ddf31SBarry Smith           dvalue++;
289e26ddf31SBarry Smith           n++;
290e26ddf31SBarry Smith           ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
291e26ddf31SBarry Smith         }
2928c74ee41SBarry Smith         ierr = PetscTokenDestroy(&token);CHKERRQ(ierr);
293e26ddf31SBarry Smith       }
294e26ddf31SBarry Smith       break;
2956356e834SBarry Smith     case OPTION_INT:
296e55864a3SBarry 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);
2973fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
2983fc1eb6aSBarry Smith       if (str[0]) {
299d25d7f95SJed Brown #if defined(PETSC_SIZEOF_LONG_LONG)
300d25d7f95SJed Brown         long long lid;
301d25d7f95SJed Brown         sscanf(str,"%lld",&lid);
3021a1499c8SBarry Smith         if (lid > PETSC_MAX_INT || lid < PETSC_MIN_INT) SETERRQ3(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Argument: -%s%s %lld",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",next->option+1,lid);
303c272547aSJed Brown #else
304d25d7f95SJed Brown         long  lid;
305d25d7f95SJed Brown         sscanf(str,"%ld",&lid);
3061a1499c8SBarry Smith         if (lid > PETSC_MAX_INT || lid < PETSC_MIN_INT) SETERRQ3(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Argument: -%s%s %ld",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",next->option+1,lid);
307c272547aSJed Brown #endif
308a297a907SKarl Rupp 
309d25d7f95SJed Brown         next->set = PETSC_TRUE;
310d25d7f95SJed Brown         *((PetscInt*)next->data) = (PetscInt)lid;
311aee2cecaSBarry Smith       }
312aee2cecaSBarry Smith       break;
313aee2cecaSBarry Smith     case OPTION_REAL:
314e55864a3SBarry 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);
3153fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3163fc1eb6aSBarry Smith       if (str[0]) {
317ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
318a4404d99SBarry Smith         sscanf(str,"%e",&ir);
319ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL_DOUBLE)
320aee2cecaSBarry Smith         sscanf(str,"%le",&ir);
321ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL___FLOAT128)
322d9822059SBarry Smith         ir = strtoflt128(str,0);
323d9822059SBarry Smith #else
324513dbe71SLisandro Dalcin         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Unknown scalar type");
325a4404d99SBarry Smith #endif
326aee2cecaSBarry Smith         next->set                 = PETSC_TRUE;
327aee2cecaSBarry Smith         *((PetscReal*)next->data) = ir;
328aee2cecaSBarry Smith       }
329aee2cecaSBarry Smith       break;
3307781c08eSBarry Smith     case OPTION_BOOL:
33183355fc5SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <%s>: %s (%s) ",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",next->option+1,*(PetscBool*)next->data ? "true": "false",next->text,next->man);CHKERRQ(ierr);
3327781c08eSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3337781c08eSBarry Smith       if (str[0]) {
3347781c08eSBarry Smith         ierr = PetscOptionsStringToBool(str,&bid);CHKERRQ(ierr);
3357781c08eSBarry Smith         next->set = PETSC_TRUE;
3367781c08eSBarry Smith         *((PetscBool*)next->data) = bid;
3377781c08eSBarry Smith       }
3387781c08eSBarry Smith       break;
339aee2cecaSBarry Smith     case OPTION_STRING:
340e55864a3SBarry 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);
3413fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3423fc1eb6aSBarry Smith       if (str[0]) {
343aee2cecaSBarry Smith         next->set = PETSC_TRUE;
34464facd6cSBarry Smith         /* must use system malloc since SAWs may free this */
3455b02f95dSBarry Smith         ierr = PetscStrdup(str,(char**)&next->data);CHKERRQ(ierr);
3466356e834SBarry Smith       }
3476356e834SBarry Smith       break;
348a264d7a6SBarry Smith     case OPTION_FLIST:
349e55864a3SBarry Smith       ierr = PetscFunctionListPrintTypes(PETSC_COMM_WORLD,stdout,PetscOptionsObject->prefix,next->option,next->text,next->man,next->flist,(char*)next->data);CHKERRQ(ierr);
3503cc1e11dSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3513cc1e11dSBarry Smith       if (str[0]) {
352e55864a3SBarry Smith         PetscOptionsObject->changedmethod = PETSC_TRUE;
3533cc1e11dSBarry Smith         next->set = PETSC_TRUE;
35464facd6cSBarry Smith         /* must use system malloc since SAWs may free this */
3555b02f95dSBarry Smith         ierr = PetscStrdup(str,(char**)&next->data);CHKERRQ(ierr);
3563cc1e11dSBarry Smith       }
3573cc1e11dSBarry Smith       break;
358b432afdaSMatthew Knepley     default:
359b432afdaSMatthew Knepley       break;
3606356e834SBarry Smith     }
3616356e834SBarry Smith     next = next->next;
3626356e834SBarry Smith   }
3636356e834SBarry Smith   PetscFunctionReturn(0);
3646356e834SBarry Smith }
3656356e834SBarry Smith 
366e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
367e04113cfSBarry Smith #include <petscviewersaws.h>
368d5649816SBarry Smith 
369d5649816SBarry Smith static int count = 0;
370d5649816SBarry Smith 
371b3506946SBarry Smith #undef __FUNCT__
372e04113cfSBarry Smith #define __FUNCT__ "PetscOptionsSAWsDestroy"
373e04113cfSBarry Smith PetscErrorCode PetscOptionsSAWsDestroy(void)
374d5649816SBarry Smith {
3752657e9d9SBarry Smith   PetscFunctionBegin;
376d5649816SBarry Smith   PetscFunctionReturn(0);
377d5649816SBarry Smith }
378d5649816SBarry Smith 
3799c1e0ce8SBarry Smith static const char *OptionsHeader = "<head>\n"
38023a1ff79SBarry Smith                                    "<script type=\"text/javascript\" src=\"http://www.mcs.anl.gov/research/projects/saws/js/jquery-1.9.1.js\"></script>\n"
38123a1ff79SBarry Smith                                    "<script type=\"text/javascript\" src=\"http://www.mcs.anl.gov/research/projects/saws/js/SAWs.js\"></script>\n"
382d1fc0251SBarry Smith                                    "<script type=\"text/javascript\" src=\"js/PETSc.js\"></script>\n"
38364bbc9efSBarry Smith                                    "<script>\n"
38464bbc9efSBarry Smith                                       "jQuery(document).ready(function() {\n"
38564bbc9efSBarry Smith                                          "PETSc.getAndDisplayDirectory(null,\"#variablesInfo\")\n"
38664bbc9efSBarry Smith                                       "})\n"
38764bbc9efSBarry Smith                                   "</script>\n"
38864bbc9efSBarry Smith                                   "</head>\n";
3891423471aSBarry Smith 
3901423471aSBarry Smith /*  Determines the size and style of the scroll region where PETSc options selectable from users are displayed */
3911423471aSBarry Smith static const char *OptionsBodyBottom = "<div id=\"variablesInfo\" style=\"background-color:lightblue;height:auto;max-height:500px;overflow:scroll;\"></div>\n<br>\n</body>";
39264bbc9efSBarry Smith 
393d5649816SBarry Smith #undef __FUNCT__
3947781c08eSBarry Smith #define __FUNCT__ "PetscOptionsSAWsInput"
395b3506946SBarry Smith /*
3967781c08eSBarry Smith     PetscOptionsSAWsInput - Presents all the PETSc Options processed by the program so the user may change them at runtime using the SAWs
397b3506946SBarry Smith 
398b3506946SBarry Smith     Bugs:
399b3506946SBarry Smith +    All processes must traverse through the exact same set of option queries do to the call to PetscScanString()
400b3506946SBarry Smith .    Internal strings have arbitrary length and string copies are not checked that they fit into string space
401b3506946SBarry Smith -    Only works for PetscInt == int, PetscReal == double etc
402b3506946SBarry Smith 
403b3506946SBarry Smith 
404b3506946SBarry Smith */
4054416b707SBarry Smith PetscErrorCode PetscOptionsSAWsInput(PetscOptionItems *PetscOptionsObject)
406b3506946SBarry Smith {
407b3506946SBarry Smith   PetscErrorCode  ierr;
4084416b707SBarry Smith   PetscOptionItem next     = PetscOptionsObject->next;
409d5649816SBarry Smith   static int      mancount = 0;
410b3506946SBarry Smith   char            options[16];
4117aab2a10SBarry Smith   PetscBool       changedmethod = PETSC_FALSE;
412a23eb982SSurtai Han   PetscBool       stopasking    = PETSC_FALSE;
41388a9752cSBarry Smith   char            manname[16],textname[16];
4142657e9d9SBarry Smith   char            dir[1024];
415b3506946SBarry Smith 
4162a409bb0SBarry Smith   PetscFunctionBegin;
417b3506946SBarry Smith   /* the next line is a bug, this will only work if all processors are here, the comm passed in is ignored!!! */
418b3506946SBarry Smith   sprintf(options,"Options_%d",count++);
419a297a907SKarl Rupp 
4207325c4a4SBarry Smith   PetscOptionsObject->pprefix = PetscOptionsObject->prefix; /* SAWs will change this, so cannot pass prefix directly */
4211bc75a8dSBarry Smith 
422d50831c4SBarry Smith   ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s","_title");CHKERRQ(ierr);
42383355fc5SBarry Smith   PetscStackCallSAWs(SAWs_Register,(dir,&PetscOptionsObject->title,1,SAWs_READ,SAWs_STRING));
4247781c08eSBarry Smith   ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s","prefix");CHKERRQ(ierr);
42583355fc5SBarry Smith   PetscStackCallSAWs(SAWs_Register,(dir,&PetscOptionsObject->pprefix,1,SAWs_READ,SAWs_STRING));
4262657e9d9SBarry Smith   PetscStackCallSAWs(SAWs_Register,("/PETSc/Options/ChangedMethod",&changedmethod,1,SAWs_WRITE,SAWs_BOOLEAN));
427a23eb982SSurtai Han   PetscStackCallSAWs(SAWs_Register,("/PETSc/Options/StopAsking",&stopasking,1,SAWs_WRITE,SAWs_BOOLEAN));
428b3506946SBarry Smith 
429b3506946SBarry Smith   while (next) {
430d50831c4SBarry Smith     sprintf(manname,"_man_%d",mancount);
4312657e9d9SBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",manname);CHKERRQ(ierr);
4327781c08eSBarry Smith     PetscStackCallSAWs(SAWs_Register,(dir,&next->man,1,SAWs_READ,SAWs_STRING));
433d50831c4SBarry Smith     sprintf(textname,"_text_%d",mancount++);
4347781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",textname);CHKERRQ(ierr);
4357781c08eSBarry Smith     PetscStackCallSAWs(SAWs_Register,(dir,&next->text,1,SAWs_READ,SAWs_STRING));
4369f32e415SBarry Smith 
437b3506946SBarry Smith     switch (next->type) {
438b3506946SBarry Smith     case OPTION_HEAD:
439b3506946SBarry Smith       break;
440b3506946SBarry Smith     case OPTION_INT_ARRAY:
4417781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4422657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_INT));
443b3506946SBarry Smith       break;
444b3506946SBarry Smith     case OPTION_REAL_ARRAY:
4457781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4462657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_DOUBLE));
447b3506946SBarry Smith       break;
448b3506946SBarry Smith     case OPTION_INT:
4497781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4502657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_INT));
451b3506946SBarry Smith       break;
452b3506946SBarry Smith     case OPTION_REAL:
4537781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4542657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_DOUBLE));
455b3506946SBarry Smith       break;
4567781c08eSBarry Smith     case OPTION_BOOL:
4577781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4582657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_BOOLEAN));
4591ae3d29cSBarry Smith       break;
4607781c08eSBarry Smith     case OPTION_BOOL_ARRAY:
4617781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4622657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_BOOLEAN));
46371f08665SBarry Smith       break;
464b3506946SBarry Smith     case OPTION_STRING:
4657781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4667781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
4671ae3d29cSBarry Smith       break;
4681ae3d29cSBarry Smith     case OPTION_STRING_ARRAY:
4697781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4702657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_STRING));
471b3506946SBarry Smith       break;
472a264d7a6SBarry Smith     case OPTION_FLIST:
473a264d7a6SBarry Smith       {
474a264d7a6SBarry Smith       PetscInt ntext;
4757781c08eSBarry Smith       ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4767781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
477a264d7a6SBarry Smith       ierr = PetscFunctionListGet(next->flist,(const char***)&next->edata,&ntext);CHKERRQ(ierr);
478a264d7a6SBarry Smith       PetscStackCallSAWs(SAWs_Set_Legal_Variable_Values,(dir,ntext,next->edata));
479a264d7a6SBarry Smith       }
480a264d7a6SBarry Smith       break;
4811ae3d29cSBarry Smith     case OPTION_ELIST:
482a264d7a6SBarry Smith       {
483a264d7a6SBarry Smith       PetscInt ntext = next->nlist;
4847781c08eSBarry Smith       ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4857781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
486ead66b60SBarry Smith       ierr = PetscMalloc1((ntext+1),(char***)&next->edata);CHKERRQ(ierr);
4871ae3d29cSBarry Smith       ierr = PetscMemcpy(next->edata,next->list,ntext*sizeof(char*));CHKERRQ(ierr);
488a264d7a6SBarry Smith       PetscStackCallSAWs(SAWs_Set_Legal_Variable_Values,(dir,ntext,next->edata));
489a264d7a6SBarry Smith       }
490a264d7a6SBarry Smith       break;
491b3506946SBarry Smith     default:
492b3506946SBarry Smith       break;
493b3506946SBarry Smith     }
494b3506946SBarry Smith     next = next->next;
495b3506946SBarry Smith   }
496b3506946SBarry Smith 
497b3506946SBarry Smith   /* wait until accessor has unlocked the memory */
49864bbc9efSBarry Smith   PetscStackCallSAWs(SAWs_Push_Header,("index.html",OptionsHeader));
49964bbc9efSBarry Smith   PetscStackCallSAWs(SAWs_Push_Body,("index.html",2,OptionsBodyBottom));
5007aab2a10SBarry Smith   ierr = PetscSAWsBlock();CHKERRQ(ierr);
50164bbc9efSBarry Smith   PetscStackCallSAWs(SAWs_Pop_Header,("index.html"));
50264bbc9efSBarry Smith   PetscStackCallSAWs(SAWs_Pop_Body,("index.html",2));
503b3506946SBarry Smith 
50488a9752cSBarry Smith   /* determine if any values have been set in GUI */
50583355fc5SBarry Smith   next = PetscOptionsObject->next;
50688a9752cSBarry Smith   while (next) {
50788a9752cSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
508f7b25cf6SBarry Smith     PetscStackCallSAWs(SAWs_Selected,(dir,(int*)&next->set));
50988a9752cSBarry Smith     next = next->next;
51088a9752cSBarry Smith   }
51188a9752cSBarry Smith 
512b3506946SBarry Smith   /* reset counter to -2; this updates the screen with the new options for the selected method */
513f7b25cf6SBarry Smith   if (changedmethod) PetscOptionsObject->count = -2;
514b3506946SBarry Smith 
515a23eb982SSurtai Han   if (stopasking) {
516a23eb982SSurtai Han     PetscOptionsPublish      = PETSC_FALSE;
51712655325SBarry Smith     PetscOptionsObject->count = 0;//do not ask for same thing again
518a23eb982SSurtai Han   }
519a23eb982SSurtai Han 
5209a492a5cSBarry Smith   PetscStackCallSAWs(SAWs_Delete,("/PETSc/Options"));
521b3506946SBarry Smith   PetscFunctionReturn(0);
522b3506946SBarry Smith }
523b3506946SBarry Smith #endif
524b3506946SBarry Smith 
5256356e834SBarry Smith #undef __FUNCT__
52653acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEnd_Private"
5274416b707SBarry Smith PetscErrorCode PetscOptionsEnd_Private(PetscOptionItems *PetscOptionsObject)
52853acd3b1SBarry Smith {
52953acd3b1SBarry Smith   PetscErrorCode  ierr;
5304416b707SBarry Smith   PetscOptionItem last;
5316356e834SBarry Smith   char            option[256],value[1024],tmp[32];
532330cf3c9SBarry Smith   size_t          j;
53353acd3b1SBarry Smith 
53453acd3b1SBarry Smith   PetscFunctionBegin;
53583355fc5SBarry Smith   if (PetscOptionsObject->next) {
53683355fc5SBarry Smith     if (!PetscOptionsObject->count) {
537a264d7a6SBarry Smith #if defined(PETSC_HAVE_SAWS)
538f7b25cf6SBarry Smith       ierr = PetscOptionsSAWsInput(PetscOptionsObject);CHKERRQ(ierr);
539b3506946SBarry Smith #else
540e55864a3SBarry Smith       ierr = PetscOptionsGetFromTextInput(PetscOptionsObject);CHKERRQ(ierr);
541b3506946SBarry Smith #endif
542aee2cecaSBarry Smith     }
543aee2cecaSBarry Smith   }
5446356e834SBarry Smith 
545e55864a3SBarry Smith   ierr = PetscFree(PetscOptionsObject->title);CHKERRQ(ierr);
5466356e834SBarry Smith 
547e26ddf31SBarry Smith   /* reset counter to -2; this updates the screen with the new options for the selected method */
548e55864a3SBarry Smith   if (PetscOptionsObject->changedmethod) PetscOptionsObject->count = -2;
5497a72a596SBarry Smith   /* reset alreadyprinted flag */
550e55864a3SBarry Smith   PetscOptionsObject->alreadyprinted = PETSC_FALSE;
551e55864a3SBarry Smith   if (PetscOptionsObject->object) PetscOptionsObject->object->optionsprinted = PETSC_TRUE;
552e55864a3SBarry Smith   PetscOptionsObject->object = NULL;
55353acd3b1SBarry Smith 
554e55864a3SBarry Smith   while (PetscOptionsObject->next) {
555e55864a3SBarry Smith     if (PetscOptionsObject->next->set) {
556e55864a3SBarry Smith       if (PetscOptionsObject->prefix) {
55753acd3b1SBarry Smith         ierr = PetscStrcpy(option,"-");CHKERRQ(ierr);
558e55864a3SBarry Smith         ierr = PetscStrcat(option,PetscOptionsObject->prefix);CHKERRQ(ierr);
559e55864a3SBarry Smith         ierr = PetscStrcat(option,PetscOptionsObject->next->option+1);CHKERRQ(ierr);
5606356e834SBarry Smith       } else {
561e55864a3SBarry Smith         ierr = PetscStrcpy(option,PetscOptionsObject->next->option);CHKERRQ(ierr);
5626356e834SBarry Smith       }
5636356e834SBarry Smith 
564e55864a3SBarry Smith       switch (PetscOptionsObject->next->type) {
5656356e834SBarry Smith       case OPTION_HEAD:
5666356e834SBarry Smith         break;
5676356e834SBarry Smith       case OPTION_INT_ARRAY:
568e55864a3SBarry Smith         sprintf(value,"%d",(int)((PetscInt*)PetscOptionsObject->next->data)[0]);
569e55864a3SBarry Smith         for (j=1; j<PetscOptionsObject->next->arraylength; j++) {
570e55864a3SBarry Smith           sprintf(tmp,"%d",(int)((PetscInt*)PetscOptionsObject->next->data)[j]);
5716356e834SBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5726356e834SBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5736356e834SBarry Smith         }
5746356e834SBarry Smith         break;
5756356e834SBarry Smith       case OPTION_INT:
576e55864a3SBarry Smith         sprintf(value,"%d",(int) *(PetscInt*)PetscOptionsObject->next->data);
5776356e834SBarry Smith         break;
5786356e834SBarry Smith       case OPTION_REAL:
579e55864a3SBarry Smith         sprintf(value,"%g",(double) *(PetscReal*)PetscOptionsObject->next->data);
5806356e834SBarry Smith         break;
5816356e834SBarry Smith       case OPTION_REAL_ARRAY:
582e55864a3SBarry Smith         sprintf(value,"%g",(double)((PetscReal*)PetscOptionsObject->next->data)[0]);
583e55864a3SBarry Smith         for (j=1; j<PetscOptionsObject->next->arraylength; j++) {
584e55864a3SBarry Smith           sprintf(tmp,"%g",(double)((PetscReal*)PetscOptionsObject->next->data)[j]);
5856356e834SBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5866356e834SBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5876356e834SBarry Smith         }
5886356e834SBarry Smith         break;
589050cccc3SHong Zhang       case OPTION_SCALAR_ARRAY:
59095f3a755SJose E. Roman         sprintf(value,"%g+%gi",(double)PetscRealPart(((PetscScalar*)PetscOptionsObject->next->data)[0]),(double)PetscImaginaryPart(((PetscScalar*)PetscOptionsObject->next->data)[0]));
591050cccc3SHong Zhang         for (j=1; j<PetscOptionsObject->next->arraylength; j++) {
59295f3a755SJose E. Roman           sprintf(tmp,"%g+%gi",(double)PetscRealPart(((PetscScalar*)PetscOptionsObject->next->data)[j]),(double)PetscImaginaryPart(((PetscScalar*)PetscOptionsObject->next->data)[j]));
593050cccc3SHong Zhang           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
594050cccc3SHong Zhang           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
595050cccc3SHong Zhang         }
596050cccc3SHong Zhang         break;
5977781c08eSBarry Smith       case OPTION_BOOL:
598e55864a3SBarry Smith         sprintf(value,"%d",*(int*)PetscOptionsObject->next->data);
5996356e834SBarry Smith         break;
6007781c08eSBarry Smith       case OPTION_BOOL_ARRAY:
601e55864a3SBarry Smith         sprintf(value,"%d",(int)((PetscBool*)PetscOptionsObject->next->data)[0]);
602e55864a3SBarry Smith         for (j=1; j<PetscOptionsObject->next->arraylength; j++) {
603e55864a3SBarry Smith           sprintf(tmp,"%d",(int)((PetscBool*)PetscOptionsObject->next->data)[j]);
6041ae3d29cSBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
6051ae3d29cSBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
6061ae3d29cSBarry Smith         }
6071ae3d29cSBarry Smith         break;
608a264d7a6SBarry Smith       case OPTION_FLIST:
6096991f827SBarry Smith         ierr = PetscStrcpy(value,(char*)PetscOptionsObject->next->data);CHKERRQ(ierr);
6106991f827SBarry Smith         break;
6111ae3d29cSBarry Smith       case OPTION_ELIST:
612e55864a3SBarry Smith         ierr = PetscStrcpy(value,(char*)PetscOptionsObject->next->data);CHKERRQ(ierr);
6136356e834SBarry Smith         break;
6141ae3d29cSBarry Smith       case OPTION_STRING:
615e55864a3SBarry Smith         ierr = PetscStrcpy(value,(char*)PetscOptionsObject->next->data);CHKERRQ(ierr);
616d51da6bfSBarry Smith         break;
6171ae3d29cSBarry Smith       case OPTION_STRING_ARRAY:
618e55864a3SBarry Smith         sprintf(value,"%s",((char**)PetscOptionsObject->next->data)[0]);
619e55864a3SBarry Smith         for (j=1; j<PetscOptionsObject->next->arraylength; j++) {
620e55864a3SBarry Smith           sprintf(tmp,"%s",((char**)PetscOptionsObject->next->data)[j]);
6211ae3d29cSBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
6221ae3d29cSBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
6231ae3d29cSBarry Smith         }
6246356e834SBarry Smith         break;
6256356e834SBarry Smith       }
626c5929fdfSBarry Smith       ierr = PetscOptionsSetValue(PetscOptionsObject->options,option,value);CHKERRQ(ierr);
6276356e834SBarry Smith     }
6286991f827SBarry Smith     if (PetscOptionsObject->next->type == OPTION_ELIST) {
6296991f827SBarry Smith       ierr = PetscStrNArrayDestroy(PetscOptionsObject->next->nlist,(char ***)&PetscOptionsObject->next->list);CHKERRQ(ierr);
6306991f827SBarry Smith     }
631e55864a3SBarry Smith     ierr   = PetscFree(PetscOptionsObject->next->text);CHKERRQ(ierr);
632e55864a3SBarry Smith     ierr   = PetscFree(PetscOptionsObject->next->option);CHKERRQ(ierr);
633e55864a3SBarry Smith     ierr   = PetscFree(PetscOptionsObject->next->man);CHKERRQ(ierr);
634e55864a3SBarry Smith     ierr   = PetscFree(PetscOptionsObject->next->edata);CHKERRQ(ierr);
635c979a496SBarry Smith 
63683355fc5SBarry Smith     if ((PetscOptionsObject->next->type == OPTION_STRING) || (PetscOptionsObject->next->type == OPTION_FLIST) || (PetscOptionsObject->next->type == OPTION_ELIST)){
63783355fc5SBarry Smith       free(PetscOptionsObject->next->data);
638c979a496SBarry Smith     } else {
63983355fc5SBarry Smith       ierr   = PetscFree(PetscOptionsObject->next->data);CHKERRQ(ierr);
640c979a496SBarry Smith     }
6417781c08eSBarry Smith 
64283355fc5SBarry Smith     last                    = PetscOptionsObject->next;
64383355fc5SBarry Smith     PetscOptionsObject->next = PetscOptionsObject->next->next;
6446356e834SBarry Smith     ierr                    = PetscFree(last);CHKERRQ(ierr);
6456356e834SBarry Smith   }
646f59f755dSBarry Smith   ierr = PetscFree(PetscOptionsObject->prefix);CHKERRQ(ierr);
647e55864a3SBarry Smith   PetscOptionsObject->next = 0;
64853acd3b1SBarry Smith   PetscFunctionReturn(0);
64953acd3b1SBarry Smith }
65053acd3b1SBarry Smith 
65153acd3b1SBarry Smith #undef __FUNCT__
652e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsEnum_Private"
65353acd3b1SBarry Smith /*@C
65453acd3b1SBarry Smith    PetscOptionsEnum - Gets the enum value for a particular option in the database.
65553acd3b1SBarry Smith 
6563f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
65753acd3b1SBarry Smith 
65853acd3b1SBarry Smith    Input Parameters:
65953acd3b1SBarry Smith +  opt - option name
66053acd3b1SBarry Smith .  text - short string that describes the option
66153acd3b1SBarry Smith .  man - manual page with additional information on option
66253acd3b1SBarry Smith .  list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null
6630fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
6640fdccdaeSBarry Smith $                 PetscOptionsEnum(..., obj->value,&object->value,...) or
6650fdccdaeSBarry Smith $                 value = defaultvalue
6660fdccdaeSBarry Smith $                 PetscOptionsEnum(..., value,&value,&flg);
6670fdccdaeSBarry Smith $                 if (flg) {
66853acd3b1SBarry Smith 
66953acd3b1SBarry Smith    Output Parameter:
67053acd3b1SBarry Smith +  value - the  value to return
671b32e0204SMatthew G Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
67253acd3b1SBarry Smith 
67353acd3b1SBarry Smith    Level: beginner
67453acd3b1SBarry Smith 
67553acd3b1SBarry Smith    Concepts: options database
67653acd3b1SBarry Smith 
67753acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
67853acd3b1SBarry Smith 
67953acd3b1SBarry Smith           list is usually something like PCASMTypes or some other predefined list of enum names
68053acd3b1SBarry Smith 
68153acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
682acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
683acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
68453acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
68553acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
686acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
687a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
68853acd3b1SBarry Smith @*/
6894416b707SBarry Smith PetscErrorCode  PetscOptionsEnum_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],const char *const *list,PetscEnum currentvalue,PetscEnum *value,PetscBool  *set)
69053acd3b1SBarry Smith {
69153acd3b1SBarry Smith   PetscErrorCode ierr;
69253acd3b1SBarry Smith   PetscInt       ntext = 0;
693aa5bb8c0SSatish Balay   PetscInt       tval;
694ace3abfcSBarry Smith   PetscBool      tflg;
69553acd3b1SBarry Smith 
69653acd3b1SBarry Smith   PetscFunctionBegin;
69753acd3b1SBarry Smith   while (list[ntext++]) {
698e32f2f54SBarry Smith     if (ntext > 50) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument appears to be wrong or have more than 50 entries");
69953acd3b1SBarry Smith   }
700e32f2f54SBarry Smith   if (ntext < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument must have at least two entries: typename and type prefix");
70153acd3b1SBarry Smith   ntext -= 3;
702e55864a3SBarry Smith   ierr   = PetscOptionsEList_Private(PetscOptionsObject,opt,text,man,list,ntext,list[currentvalue],&tval,&tflg);CHKERRQ(ierr);
703aa5bb8c0SSatish Balay   /* with PETSC_USE_64BIT_INDICES sizeof(PetscInt) != sizeof(PetscEnum) */
704aa5bb8c0SSatish Balay   if (tflg) *value = (PetscEnum)tval;
705aa5bb8c0SSatish Balay   if (set)  *set   = tflg;
70653acd3b1SBarry Smith   PetscFunctionReturn(0);
70753acd3b1SBarry Smith }
70853acd3b1SBarry Smith 
709d3e47460SLisandro Dalcin #undef __FUNCT__
710d3e47460SLisandro Dalcin #define __FUNCT__ "PetscOptionsEnumArray_Private"
711d3e47460SLisandro Dalcin /*@C
712d3e47460SLisandro Dalcin    PetscOptionsEnumArray - Gets an array of enum values for a particular
713d3e47460SLisandro Dalcin    option in the database.
714d3e47460SLisandro Dalcin 
715d3e47460SLisandro Dalcin    Logically Collective on the communicator passed in PetscOptionsBegin()
716d3e47460SLisandro Dalcin 
717d3e47460SLisandro Dalcin    Input Parameters:
718d3e47460SLisandro Dalcin +  opt - the option one is seeking
719d3e47460SLisandro Dalcin .  text - short string describing option
720d3e47460SLisandro Dalcin .  man - manual page for option
721d3e47460SLisandro Dalcin -  n - maximum number of values
722d3e47460SLisandro Dalcin 
723d3e47460SLisandro Dalcin    Output Parameter:
724d3e47460SLisandro Dalcin +  value - location to copy values
725d3e47460SLisandro Dalcin .  n - actual number of values found
726d3e47460SLisandro Dalcin -  set - PETSC_TRUE if found, else PETSC_FALSE
727d3e47460SLisandro Dalcin 
728d3e47460SLisandro Dalcin    Level: beginner
729d3e47460SLisandro Dalcin 
730d3e47460SLisandro Dalcin    Notes:
731d3e47460SLisandro Dalcin    The array must be passed as a comma separated list.
732d3e47460SLisandro Dalcin 
733d3e47460SLisandro Dalcin    There must be no intervening spaces between the values.
734d3e47460SLisandro Dalcin 
735d3e47460SLisandro Dalcin    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
736d3e47460SLisandro Dalcin 
737d3e47460SLisandro Dalcin    Concepts: options database^array of enums
738d3e47460SLisandro Dalcin 
739d3e47460SLisandro Dalcin .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
740d3e47460SLisandro Dalcin           PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
741d3e47460SLisandro Dalcin           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
742d3e47460SLisandro Dalcin           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
743d3e47460SLisandro Dalcin           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
744d3e47460SLisandro Dalcin           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsRealArray()
745d3e47460SLisandro Dalcin @*/
7464416b707SBarry Smith PetscErrorCode  PetscOptionsEnumArray_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],const char *const *list,PetscEnum value[],PetscInt *n,PetscBool  *set)
747d3e47460SLisandro Dalcin {
748d3e47460SLisandro Dalcin   PetscInt        i,nlist = 0;
7494416b707SBarry Smith   PetscOptionItem amsopt;
750d3e47460SLisandro Dalcin   PetscErrorCode  ierr;
751d3e47460SLisandro Dalcin 
752d3e47460SLisandro Dalcin   PetscFunctionBegin;
753d3e47460SLisandro Dalcin   while (list[nlist++]) if (nlist > 50) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument appears to be wrong or have more than 50 entries");
754d3e47460SLisandro Dalcin   if (nlist < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument must have at least two entries: typename and type prefix");
755d3e47460SLisandro Dalcin   nlist -= 3; /* drop enum name, prefix, and null termination */
756d3e47460SLisandro Dalcin   if (0 && !PetscOptionsObject->count) { /* XXX Requires additional support */
757d3e47460SLisandro Dalcin     PetscEnum *vals;
7584416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_INT_ARRAY/*XXX OPTION_ENUM_ARRAY*/,&amsopt);CHKERRQ(ierr);
759d3e47460SLisandro Dalcin     ierr = PetscStrNArrayallocpy(nlist,list,(char***)&amsopt->list);CHKERRQ(ierr);
760d3e47460SLisandro Dalcin     amsopt->nlist = nlist;
761d3e47460SLisandro Dalcin     ierr = PetscMalloc1(*n,(PetscEnum**)&amsopt->data);CHKERRQ(ierr);
762d3e47460SLisandro Dalcin     amsopt->arraylength = *n;
763d3e47460SLisandro Dalcin     vals = (PetscEnum*)amsopt->data;
764d3e47460SLisandro Dalcin     for (i=0; i<*n; i++) vals[i] = value[i];
765d3e47460SLisandro Dalcin   }
766c5929fdfSBarry Smith   ierr = PetscOptionsGetEnumArray(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,list,value,n,set);CHKERRQ(ierr);
767d3e47460SLisandro Dalcin   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
768d3e47460SLisandro Dalcin     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%s",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,list[value[0]]);CHKERRQ(ierr);
769d3e47460SLisandro Dalcin     for (i=1; i<*n; i++) {ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%s",list[value[i]]);CHKERRQ(ierr);}
770d3e47460SLisandro Dalcin     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (choose from)",text);CHKERRQ(ierr);
771d3e47460SLisandro Dalcin     for (i=0; i<nlist; i++) {ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," %s",list[i]);CHKERRQ(ierr);}
772d3e47460SLisandro Dalcin     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," (%s)\n",ManSection(man));CHKERRQ(ierr);
773d3e47460SLisandro Dalcin   }
774d3e47460SLisandro Dalcin   PetscFunctionReturn(0);
775d3e47460SLisandro Dalcin }
776d3e47460SLisandro Dalcin 
77753acd3b1SBarry Smith /* -------------------------------------------------------------------------------------------------------------*/
77853acd3b1SBarry Smith #undef __FUNCT__
779e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsInt_Private"
78053acd3b1SBarry Smith /*@C
78153acd3b1SBarry Smith    PetscOptionsInt - Gets the integer value for a particular option in the database.
78253acd3b1SBarry Smith 
7833f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
78453acd3b1SBarry Smith 
78553acd3b1SBarry Smith    Input Parameters:
78653acd3b1SBarry Smith +  opt - option name
78753acd3b1SBarry Smith .  text - short string that describes the option
78853acd3b1SBarry Smith .  man - manual page with additional information on option
7890fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
7900fdccdaeSBarry Smith $                 PetscOptionsInt(..., obj->value,&object->value,...) or
7910fdccdaeSBarry Smith $                 value = defaultvalue
7920fdccdaeSBarry Smith $                 PetscOptionsInt(..., value,&value,&flg);
7930fdccdaeSBarry Smith $                 if (flg) {
79453acd3b1SBarry Smith 
79553acd3b1SBarry Smith    Output Parameter:
79653acd3b1SBarry Smith +  value - the integer value to return
79753acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
79853acd3b1SBarry Smith 
79953acd3b1SBarry Smith    Level: beginner
80053acd3b1SBarry Smith 
80153acd3b1SBarry Smith    Concepts: options database^has int
80253acd3b1SBarry Smith 
80353acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
80453acd3b1SBarry Smith 
80553acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
806acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
807acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
80853acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
80953acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
810acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
811a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
81253acd3b1SBarry Smith @*/
8134416b707SBarry Smith PetscErrorCode  PetscOptionsInt_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscInt currentvalue,PetscInt *value,PetscBool  *set)
81453acd3b1SBarry Smith {
81553acd3b1SBarry Smith   PetscErrorCode  ierr;
8164416b707SBarry Smith   PetscOptionItem amsopt;
81712655325SBarry Smith   PetscBool       wasset;
81853acd3b1SBarry Smith 
81953acd3b1SBarry Smith   PetscFunctionBegin;
820e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
8214416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_INT,&amsopt);CHKERRQ(ierr);
8226356e834SBarry Smith     ierr = PetscMalloc(sizeof(PetscInt),&amsopt->data);CHKERRQ(ierr);
82312655325SBarry Smith     *(PetscInt*)amsopt->data = currentvalue;
8243e211508SBarry Smith 
825c5929fdfSBarry Smith     ierr = PetscOptionsGetInt(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,&currentvalue,&wasset);CHKERRQ(ierr);
8263e211508SBarry Smith     if (wasset) {
82712655325SBarry Smith       *(PetscInt*)amsopt->data = currentvalue;
8283e211508SBarry Smith     }
829af6d86caSBarry Smith   }
830c5929fdfSBarry Smith   ierr = PetscOptionsGetInt(NULL,PetscOptionsObject->prefix,opt,value,set);CHKERRQ(ierr);
831e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
8321a1499c8SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%d>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,currentvalue,text,ManSection(man));CHKERRQ(ierr);
83353acd3b1SBarry Smith   }
83453acd3b1SBarry Smith   PetscFunctionReturn(0);
83553acd3b1SBarry Smith }
83653acd3b1SBarry Smith 
83753acd3b1SBarry Smith #undef __FUNCT__
838e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsString_Private"
83953acd3b1SBarry Smith /*@C
84053acd3b1SBarry Smith    PetscOptionsString - Gets the string value for a particular option in the database.
84153acd3b1SBarry Smith 
8423f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
84353acd3b1SBarry Smith 
84453acd3b1SBarry Smith    Input Parameters:
84553acd3b1SBarry Smith +  opt - option name
84653acd3b1SBarry Smith .  text - short string that describes the option
84753acd3b1SBarry Smith .  man - manual page with additional information on option
8480fdccdaeSBarry Smith .  currentvalue - the current value; caller is responsible for setting this value correctly. This is not used to set value
849bcbf2dc5SJed Brown -  len - length of the result string including null terminator
85053acd3b1SBarry Smith 
85153acd3b1SBarry Smith    Output Parameter:
85253acd3b1SBarry Smith +  value - the value to return
85353acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
85453acd3b1SBarry Smith 
85553acd3b1SBarry Smith    Level: beginner
85653acd3b1SBarry Smith 
85753acd3b1SBarry Smith    Concepts: options database^has int
85853acd3b1SBarry Smith 
85953acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
86053acd3b1SBarry Smith 
8617fccdfe4SBarry 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).
8627fccdfe4SBarry Smith 
863c5929fdfSBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(NULL,),
864acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
865acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsReal(), PetscOptionsBool(),
86653acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
86753acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
868acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
869a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
87053acd3b1SBarry Smith @*/
8714416b707SBarry Smith PetscErrorCode  PetscOptionsString_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],const char currentvalue[],char value[],size_t len,PetscBool  *set)
87253acd3b1SBarry Smith {
87353acd3b1SBarry Smith   PetscErrorCode  ierr;
8744416b707SBarry Smith   PetscOptionItem amsopt;
87553acd3b1SBarry Smith 
87653acd3b1SBarry Smith   PetscFunctionBegin;
8771a1499c8SBarry Smith   if (!PetscOptionsObject->count) {
8784416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
87964facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
8800fdccdaeSBarry Smith     ierr = PetscStrdup(currentvalue ? currentvalue : "",(char**)&amsopt->data);CHKERRQ(ierr);
881af6d86caSBarry Smith   }
882c5929fdfSBarry Smith   ierr = PetscOptionsGetString(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,len,set);CHKERRQ(ierr);
883e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
8841a1499c8SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%s>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,currentvalue,text,ManSection(man));CHKERRQ(ierr);
88553acd3b1SBarry Smith   }
88653acd3b1SBarry Smith   PetscFunctionReturn(0);
88753acd3b1SBarry Smith }
88853acd3b1SBarry Smith 
88953acd3b1SBarry Smith #undef __FUNCT__
890e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsReal_Private"
89153acd3b1SBarry Smith /*@C
89253acd3b1SBarry Smith    PetscOptionsReal - Gets the PetscReal value for a particular option in the database.
89353acd3b1SBarry Smith 
8943f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
89553acd3b1SBarry Smith 
89653acd3b1SBarry Smith    Input Parameters:
89753acd3b1SBarry Smith +  opt - option name
89853acd3b1SBarry Smith .  text - short string that describes the option
89953acd3b1SBarry Smith .  man - manual page with additional information on option
9000fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
9010fdccdaeSBarry Smith $                 PetscOptionsReal(..., obj->value,&object->value,...) or
9020fdccdaeSBarry Smith $                 value = defaultvalue
9030fdccdaeSBarry Smith $                 PetscOptionsReal(..., value,&value,&flg);
9040fdccdaeSBarry Smith $                 if (flg) {
90553acd3b1SBarry Smith 
90653acd3b1SBarry Smith    Output Parameter:
90753acd3b1SBarry Smith +  value - the value to return
90853acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
90953acd3b1SBarry Smith 
91053acd3b1SBarry Smith    Level: beginner
91153acd3b1SBarry Smith 
91253acd3b1SBarry Smith    Concepts: options database^has int
91353acd3b1SBarry Smith 
91453acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
91553acd3b1SBarry Smith 
916c5929fdfSBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(NULL,),
917acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
918acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
91953acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
92053acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
921acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
922a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
92353acd3b1SBarry Smith @*/
9244416b707SBarry Smith PetscErrorCode  PetscOptionsReal_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscReal currentvalue,PetscReal *value,PetscBool  *set)
92553acd3b1SBarry Smith {
92653acd3b1SBarry Smith   PetscErrorCode  ierr;
9274416b707SBarry Smith   PetscOptionItem amsopt;
92853acd3b1SBarry Smith 
92953acd3b1SBarry Smith   PetscFunctionBegin;
930e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
9314416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_REAL,&amsopt);CHKERRQ(ierr);
932538aa990SBarry Smith     ierr = PetscMalloc(sizeof(PetscReal),&amsopt->data);CHKERRQ(ierr);
933a297a907SKarl Rupp 
9340fdccdaeSBarry Smith     *(PetscReal*)amsopt->data = currentvalue;
935538aa990SBarry Smith   }
936c5929fdfSBarry Smith   ierr = PetscOptionsGetReal(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,set);CHKERRQ(ierr);
9371a1499c8SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
9381a1499c8SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%g>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,(double)currentvalue,text,ManSection(man));CHKERRQ(ierr);
93953acd3b1SBarry Smith   }
94053acd3b1SBarry Smith   PetscFunctionReturn(0);
94153acd3b1SBarry Smith }
94253acd3b1SBarry Smith 
94353acd3b1SBarry Smith #undef __FUNCT__
944e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsScalar_Private"
94553acd3b1SBarry Smith /*@C
94653acd3b1SBarry Smith    PetscOptionsScalar - Gets the scalar value for a particular option in the database.
94753acd3b1SBarry Smith 
9483f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
94953acd3b1SBarry Smith 
95053acd3b1SBarry Smith    Input Parameters:
95153acd3b1SBarry Smith +  opt - option name
95253acd3b1SBarry Smith .  text - short string that describes the option
95353acd3b1SBarry Smith .  man - manual page with additional information on option
9540fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
9550fdccdaeSBarry Smith $                 PetscOptionsScalar(..., obj->value,&object->value,...) or
9560fdccdaeSBarry Smith $                 value = defaultvalue
9570fdccdaeSBarry Smith $                 PetscOptionsScalar(..., value,&value,&flg);
9580fdccdaeSBarry Smith $                 if (flg) {
9590fdccdaeSBarry Smith 
96053acd3b1SBarry Smith 
96153acd3b1SBarry Smith    Output Parameter:
96253acd3b1SBarry Smith +  value - the value to return
96353acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
96453acd3b1SBarry Smith 
96553acd3b1SBarry Smith    Level: beginner
96653acd3b1SBarry Smith 
96753acd3b1SBarry Smith    Concepts: options database^has int
96853acd3b1SBarry Smith 
96953acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
97053acd3b1SBarry Smith 
971c5929fdfSBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(NULL,),
972acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
973acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
97453acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
97553acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
976acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
977a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
97853acd3b1SBarry Smith @*/
9794416b707SBarry Smith PetscErrorCode  PetscOptionsScalar_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscScalar currentvalue,PetscScalar *value,PetscBool  *set)
98053acd3b1SBarry Smith {
98153acd3b1SBarry Smith   PetscErrorCode ierr;
98253acd3b1SBarry Smith 
98353acd3b1SBarry Smith   PetscFunctionBegin;
98453acd3b1SBarry Smith #if !defined(PETSC_USE_COMPLEX)
9850fdccdaeSBarry Smith   ierr = PetscOptionsReal(opt,text,man,currentvalue,value,set);CHKERRQ(ierr);
98653acd3b1SBarry Smith #else
987c5929fdfSBarry Smith   ierr = PetscOptionsGetScalar(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,set);CHKERRQ(ierr);
98853acd3b1SBarry Smith #endif
98953acd3b1SBarry Smith   PetscFunctionReturn(0);
99053acd3b1SBarry Smith }
99153acd3b1SBarry Smith 
99253acd3b1SBarry Smith #undef __FUNCT__
993e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsName_Private"
99453acd3b1SBarry Smith /*@C
99590d69ab7SBarry 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
99690d69ab7SBarry Smith                       its value is set to false.
99753acd3b1SBarry Smith 
9983f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
99953acd3b1SBarry Smith 
100053acd3b1SBarry Smith    Input Parameters:
100153acd3b1SBarry Smith +  opt - option name
100253acd3b1SBarry Smith .  text - short string that describes the option
100353acd3b1SBarry Smith -  man - manual page with additional information on option
100453acd3b1SBarry Smith 
100553acd3b1SBarry Smith    Output Parameter:
100653acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
100753acd3b1SBarry Smith 
100853acd3b1SBarry Smith    Level: beginner
100953acd3b1SBarry Smith 
101053acd3b1SBarry Smith    Concepts: options database^has int
101153acd3b1SBarry Smith 
101253acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
101353acd3b1SBarry Smith 
1014c5929fdfSBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(NULL,),
1015acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
1016acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
101753acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
101853acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1019acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1020a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
102153acd3b1SBarry Smith @*/
10224416b707SBarry Smith PetscErrorCode  PetscOptionsName_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool  *flg)
102353acd3b1SBarry Smith {
102453acd3b1SBarry Smith   PetscErrorCode  ierr;
10254416b707SBarry Smith   PetscOptionItem amsopt;
102653acd3b1SBarry Smith 
102753acd3b1SBarry Smith   PetscFunctionBegin;
1028e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
10294416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1030ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1031a297a907SKarl Rupp 
1032ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
10331ae3d29cSBarry Smith   }
1034c5929fdfSBarry Smith   ierr = PetscOptionsHasName(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,flg);CHKERRQ(ierr);
1035e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1036e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
103753acd3b1SBarry Smith   }
103853acd3b1SBarry Smith   PetscFunctionReturn(0);
103953acd3b1SBarry Smith }
104053acd3b1SBarry Smith 
104153acd3b1SBarry Smith #undef __FUNCT__
1042e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsFList_Private"
104353acd3b1SBarry Smith /*@C
1044a264d7a6SBarry Smith      PetscOptionsFList - Puts a list of option values that a single one may be selected from
104553acd3b1SBarry Smith 
10463f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
104753acd3b1SBarry Smith 
104853acd3b1SBarry Smith    Input Parameters:
104953acd3b1SBarry Smith +  opt - option name
105053acd3b1SBarry Smith .  text - short string that describes the option
105153acd3b1SBarry Smith .  man - manual page with additional information on option
105253acd3b1SBarry Smith .  list - the possible choices
10530fdccdaeSBarry Smith .  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with
10540fdccdaeSBarry Smith $                 PetscOptionsFlist(..., obj->value,value,len,&flg);
10550fdccdaeSBarry Smith $                 if (flg) {
10563cc1e11dSBarry Smith -  len - the length of the character array value
105753acd3b1SBarry Smith 
105853acd3b1SBarry Smith    Output Parameter:
105953acd3b1SBarry Smith +  value - the value to return
106053acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
106153acd3b1SBarry Smith 
106253acd3b1SBarry Smith    Level: intermediate
106353acd3b1SBarry Smith 
106453acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
106553acd3b1SBarry Smith 
106653acd3b1SBarry Smith    See PetscOptionsEList() for when the choices are given in a string array
106753acd3b1SBarry Smith 
106853acd3b1SBarry Smith    To get a listing of all currently specified options,
106988c29154SBarry Smith     see PetscOptionsView() or PetscOptionsGetAll()
107053acd3b1SBarry Smith 
1071eabe10d7SBarry Smith    Developer Note: This cannot check for invalid selection because of things like MATAIJ that are not included in the list
1072eabe10d7SBarry Smith 
107353acd3b1SBarry Smith    Concepts: options database^list
107453acd3b1SBarry Smith 
1075c5929fdfSBarry Smith .seealso: PetscOptionsGetInt(NULL,), PetscOptionsGetReal(),
1076acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
107753acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
107853acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1079acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1080a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsEnum()
108153acd3b1SBarry Smith @*/
10824416b707SBarry Smith PetscErrorCode  PetscOptionsFList_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char ltext[],const char man[],PetscFunctionList list,const char currentvalue[],char value[],size_t len,PetscBool  *set)
108353acd3b1SBarry Smith {
108453acd3b1SBarry Smith   PetscErrorCode  ierr;
10854416b707SBarry Smith   PetscOptionItem amsopt;
108653acd3b1SBarry Smith 
108753acd3b1SBarry Smith   PetscFunctionBegin;
10881a1499c8SBarry Smith   if (!PetscOptionsObject->count) {
10894416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,ltext,man,OPTION_FLIST,&amsopt);CHKERRQ(ierr);
109064facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
10910fdccdaeSBarry Smith     ierr = PetscStrdup(currentvalue ? currentvalue : "",(char**)&amsopt->data);CHKERRQ(ierr);
10923cc1e11dSBarry Smith     amsopt->flist = list;
10933cc1e11dSBarry Smith   }
1094c5929fdfSBarry Smith   ierr = PetscOptionsGetString(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,len,set);CHKERRQ(ierr);
10951a1499c8SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
10961a1499c8SBarry Smith     ierr = PetscFunctionListPrintTypes(PetscOptionsObject->comm,stdout,PetscOptionsObject->prefix,opt,ltext,man,list,currentvalue);CHKERRQ(ierr);CHKERRQ(ierr);
109753acd3b1SBarry Smith   }
109853acd3b1SBarry Smith   PetscFunctionReturn(0);
109953acd3b1SBarry Smith }
110053acd3b1SBarry Smith 
110153acd3b1SBarry Smith #undef __FUNCT__
1102e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsEList_Private"
110353acd3b1SBarry Smith /*@C
110453acd3b1SBarry Smith      PetscOptionsEList - Puts a list of option values that a single one may be selected from
110553acd3b1SBarry Smith 
11063f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
110753acd3b1SBarry Smith 
110853acd3b1SBarry Smith    Input Parameters:
110953acd3b1SBarry Smith +  opt - option name
111053acd3b1SBarry Smith .  ltext - short string that describes the option
111153acd3b1SBarry Smith .  man - manual page with additional information on option
1112a264d7a6SBarry Smith .  list - the possible choices (one of these must be selected, anything else is invalid)
111353acd3b1SBarry Smith .  ntext - number of choices
11140fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with
11150fdccdaeSBarry Smith $                 PetscOptionsElist(..., obj->value,&value,&flg);
11160fdccdaeSBarry Smith $                 if (flg) {
11170fdccdaeSBarry Smith 
111853acd3b1SBarry Smith 
111953acd3b1SBarry Smith    Output Parameter:
112053acd3b1SBarry Smith +  value - the index of the value to return
112153acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
112253acd3b1SBarry Smith 
112353acd3b1SBarry Smith    Level: intermediate
112453acd3b1SBarry Smith 
112553acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
112653acd3b1SBarry Smith 
1127a264d7a6SBarry Smith    See PetscOptionsFList() for when the choices are given in a PetscFunctionList()
112853acd3b1SBarry Smith 
112953acd3b1SBarry Smith    Concepts: options database^list
113053acd3b1SBarry Smith 
1131c5929fdfSBarry Smith .seealso: PetscOptionsGetInt(NULL,), PetscOptionsGetReal(),
1132acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
113353acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
113453acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1135acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1136a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEnum()
113753acd3b1SBarry Smith @*/
11384416b707SBarry Smith PetscErrorCode  PetscOptionsEList_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char ltext[],const char man[],const char *const *list,PetscInt ntext,const char currentvalue[],PetscInt *value,PetscBool  *set)
113953acd3b1SBarry Smith {
114053acd3b1SBarry Smith   PetscErrorCode  ierr;
114153acd3b1SBarry Smith   PetscInt        i;
11424416b707SBarry Smith   PetscOptionItem amsopt;
114353acd3b1SBarry Smith 
114453acd3b1SBarry Smith   PetscFunctionBegin;
11451a1499c8SBarry Smith   if (!PetscOptionsObject->count) {
11464416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,ltext,man,OPTION_ELIST,&amsopt);CHKERRQ(ierr);
114764facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
11480fdccdaeSBarry Smith     ierr = PetscStrdup(currentvalue ? currentvalue : "",(char**)&amsopt->data);CHKERRQ(ierr);
11496991f827SBarry Smith     ierr = PetscStrNArrayallocpy(ntext,list,(char***)&amsopt->list);CHKERRQ(ierr);
11501ae3d29cSBarry Smith     amsopt->nlist = ntext;
11511ae3d29cSBarry Smith   }
1152c5929fdfSBarry Smith   ierr = PetscOptionsGetEList(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,list,ntext,value,set);CHKERRQ(ierr);
11531a1499c8SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
11541a1499c8SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%s> (choose one of)",PetscOptionsObject->prefix?PetscOptionsObject->prefix:"",opt+1,currentvalue);CHKERRQ(ierr);
115553acd3b1SBarry Smith     for (i=0; i<ntext; i++) {
1156e55864a3SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," %s",list[i]);CHKERRQ(ierr);
115753acd3b1SBarry Smith     }
1158e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," (%s)\n",ManSection(man));CHKERRQ(ierr);
115953acd3b1SBarry Smith   }
116053acd3b1SBarry Smith   PetscFunctionReturn(0);
116153acd3b1SBarry Smith }
116253acd3b1SBarry Smith 
116353acd3b1SBarry Smith #undef __FUNCT__
1164e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsBoolGroupBegin_Private"
116553acd3b1SBarry Smith /*@C
1166acfcf0e5SJed Brown      PetscOptionsBoolGroupBegin - First in a series of logical queries on the options database for
1167d5649816SBarry Smith        which at most a single value can be true.
116853acd3b1SBarry Smith 
11693f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
117053acd3b1SBarry Smith 
117153acd3b1SBarry Smith    Input Parameters:
117253acd3b1SBarry Smith +  opt - option name
117353acd3b1SBarry Smith .  text - short string that describes the option
117453acd3b1SBarry Smith -  man - manual page with additional information on option
117553acd3b1SBarry Smith 
117653acd3b1SBarry Smith    Output Parameter:
117753acd3b1SBarry Smith .  flg - whether that option was set or not
117853acd3b1SBarry Smith 
117953acd3b1SBarry Smith    Level: intermediate
118053acd3b1SBarry Smith 
118153acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
118253acd3b1SBarry Smith 
1183acfcf0e5SJed Brown    Must be followed by 0 or more PetscOptionsBoolGroup()s and PetscOptionsBoolGroupEnd()
118453acd3b1SBarry Smith 
118553acd3b1SBarry Smith     Concepts: options database^logical group
118653acd3b1SBarry Smith 
1187c5929fdfSBarry Smith .seealso: PetscOptionsGetInt(NULL,), PetscOptionsGetReal(),
1188acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
118953acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
119053acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1191acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1192a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
119353acd3b1SBarry Smith @*/
11944416b707SBarry Smith PetscErrorCode  PetscOptionsBoolGroupBegin_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool  *flg)
119553acd3b1SBarry Smith {
119653acd3b1SBarry Smith   PetscErrorCode  ierr;
11974416b707SBarry Smith   PetscOptionItem amsopt;
119853acd3b1SBarry Smith 
119953acd3b1SBarry Smith   PetscFunctionBegin;
1200e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
12014416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1202ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1203a297a907SKarl Rupp 
1204ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
12051ae3d29cSBarry Smith   }
120668b16fdaSBarry Smith   *flg = PETSC_FALSE;
1207c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,flg,NULL);CHKERRQ(ierr);
1208e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1209e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  Pick at most one of -------------\n");CHKERRQ(ierr);
1210e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"    -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
121153acd3b1SBarry Smith   }
121253acd3b1SBarry Smith   PetscFunctionReturn(0);
121353acd3b1SBarry Smith }
121453acd3b1SBarry Smith 
121553acd3b1SBarry Smith #undef __FUNCT__
1216e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsBoolGroup_Private"
121753acd3b1SBarry Smith /*@C
1218acfcf0e5SJed Brown      PetscOptionsBoolGroup - One in a series of logical queries on the options database for
1219d5649816SBarry Smith        which at most a single value can be true.
122053acd3b1SBarry Smith 
12213f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
122253acd3b1SBarry Smith 
122353acd3b1SBarry Smith    Input Parameters:
122453acd3b1SBarry Smith +  opt - option name
122553acd3b1SBarry Smith .  text - short string that describes the option
122653acd3b1SBarry Smith -  man - manual page with additional information on option
122753acd3b1SBarry Smith 
122853acd3b1SBarry Smith    Output Parameter:
122953acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
123053acd3b1SBarry Smith 
123153acd3b1SBarry Smith    Level: intermediate
123253acd3b1SBarry Smith 
123353acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
123453acd3b1SBarry Smith 
1235acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin() and preceded a PetscOptionsBoolGroupEnd()
123653acd3b1SBarry Smith 
123753acd3b1SBarry Smith     Concepts: options database^logical group
123853acd3b1SBarry Smith 
1239c5929fdfSBarry Smith .seealso: PetscOptionsGetInt(NULL,), PetscOptionsGetReal(),
1240acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
124153acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
124253acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1243acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1244a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
124553acd3b1SBarry Smith @*/
12464416b707SBarry Smith PetscErrorCode  PetscOptionsBoolGroup_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool  *flg)
124753acd3b1SBarry Smith {
124853acd3b1SBarry Smith   PetscErrorCode  ierr;
12494416b707SBarry Smith   PetscOptionItem amsopt;
125053acd3b1SBarry Smith 
125153acd3b1SBarry Smith   PetscFunctionBegin;
1252e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
12534416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1254ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1255a297a907SKarl Rupp 
1256ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
12571ae3d29cSBarry Smith   }
125817326d04SJed Brown   *flg = PETSC_FALSE;
1259c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,flg,NULL);CHKERRQ(ierr);
1260e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1261e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"    -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
126253acd3b1SBarry Smith   }
126353acd3b1SBarry Smith   PetscFunctionReturn(0);
126453acd3b1SBarry Smith }
126553acd3b1SBarry Smith 
126653acd3b1SBarry Smith #undef __FUNCT__
1267e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsBoolGroupEnd_Private"
126853acd3b1SBarry Smith /*@C
1269acfcf0e5SJed Brown      PetscOptionsBoolGroupEnd - Last in a series of logical queries on the options database for
1270d5649816SBarry Smith        which at most a single value can be true.
127153acd3b1SBarry Smith 
12723f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
127353acd3b1SBarry Smith 
127453acd3b1SBarry Smith    Input Parameters:
127553acd3b1SBarry Smith +  opt - option name
127653acd3b1SBarry Smith .  text - short string that describes the option
127753acd3b1SBarry Smith -  man - manual page with additional information on option
127853acd3b1SBarry Smith 
127953acd3b1SBarry Smith    Output Parameter:
128053acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
128153acd3b1SBarry Smith 
128253acd3b1SBarry Smith    Level: intermediate
128353acd3b1SBarry Smith 
128453acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
128553acd3b1SBarry Smith 
1286acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin()
128753acd3b1SBarry Smith 
128853acd3b1SBarry Smith     Concepts: options database^logical group
128953acd3b1SBarry Smith 
1290c5929fdfSBarry Smith .seealso: PetscOptionsGetInt(NULL,), PetscOptionsGetReal(),
1291acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
129253acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
129353acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1294acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1295a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
129653acd3b1SBarry Smith @*/
12974416b707SBarry Smith PetscErrorCode  PetscOptionsBoolGroupEnd_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool  *flg)
129853acd3b1SBarry Smith {
129953acd3b1SBarry Smith   PetscErrorCode  ierr;
13004416b707SBarry Smith   PetscOptionItem amsopt;
130153acd3b1SBarry Smith 
130253acd3b1SBarry Smith   PetscFunctionBegin;
1303e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
13044416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1305ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1306a297a907SKarl Rupp 
1307ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
13081ae3d29cSBarry Smith   }
130917326d04SJed Brown   *flg = PETSC_FALSE;
1310c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,flg,NULL);CHKERRQ(ierr);
1311e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1312e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"    -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
131353acd3b1SBarry Smith   }
131453acd3b1SBarry Smith   PetscFunctionReturn(0);
131553acd3b1SBarry Smith }
131653acd3b1SBarry Smith 
131753acd3b1SBarry Smith #undef __FUNCT__
1318e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsBool_Private"
131953acd3b1SBarry Smith /*@C
1320acfcf0e5SJed Brown    PetscOptionsBool - Determines if a particular option is in the database with a true or false
132153acd3b1SBarry Smith 
13223f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
132353acd3b1SBarry Smith 
132453acd3b1SBarry Smith    Input Parameters:
132553acd3b1SBarry Smith +  opt - option name
132653acd3b1SBarry Smith .  text - short string that describes the option
1327868c398cSBarry Smith .  man - manual page with additional information on option
132894ae4db5SBarry Smith -  currentvalue - the current value
132953acd3b1SBarry Smith 
133053acd3b1SBarry Smith    Output Parameter:
133153acd3b1SBarry Smith .  flg - PETSC_TRUE or PETSC_FALSE
133253acd3b1SBarry Smith .  set - PETSC_TRUE if found, else PETSC_FALSE
133353acd3b1SBarry Smith 
133453acd3b1SBarry Smith    Level: beginner
133553acd3b1SBarry Smith 
133653acd3b1SBarry Smith    Concepts: options database^logical
133753acd3b1SBarry Smith 
133853acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
133953acd3b1SBarry Smith 
1340c5929fdfSBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(NULL,),
1341acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
1342acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
134353acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
134453acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1345acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1346a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
134753acd3b1SBarry Smith @*/
13484416b707SBarry Smith PetscErrorCode  PetscOptionsBool_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool currentvalue,PetscBool  *flg,PetscBool  *set)
134953acd3b1SBarry Smith {
135053acd3b1SBarry Smith   PetscErrorCode  ierr;
1351ace3abfcSBarry Smith   PetscBool       iset;
13524416b707SBarry Smith   PetscOptionItem amsopt;
135353acd3b1SBarry Smith 
135453acd3b1SBarry Smith   PetscFunctionBegin;
1355e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
13564416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1357ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1358a297a907SKarl Rupp 
135994ae4db5SBarry Smith     *(PetscBool*)amsopt->data = currentvalue;
1360af6d86caSBarry Smith   }
1361c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,flg,&iset);CHKERRQ(ierr);
136253acd3b1SBarry Smith   if (set) *set = iset;
13631a1499c8SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
136494ae4db5SBarry Smith     const char *v = PetscBools[currentvalue];
13651a1499c8SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s: <%s> %s (%s)\n",PetscOptionsObject->prefix?PetscOptionsObject->prefix:"",opt+1,v,text,ManSection(man));CHKERRQ(ierr);
136653acd3b1SBarry Smith   }
136753acd3b1SBarry Smith   PetscFunctionReturn(0);
136853acd3b1SBarry Smith }
136953acd3b1SBarry Smith 
137053acd3b1SBarry Smith #undef __FUNCT__
1371e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsRealArray_Private"
137253acd3b1SBarry Smith /*@C
137353acd3b1SBarry Smith    PetscOptionsRealArray - Gets an array of double values for a particular
137453acd3b1SBarry Smith    option in the database. The values must be separated with commas with
137553acd3b1SBarry Smith    no intervening spaces.
137653acd3b1SBarry Smith 
13773f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
137853acd3b1SBarry Smith 
137953acd3b1SBarry Smith    Input Parameters:
138053acd3b1SBarry Smith +  opt - the option one is seeking
138153acd3b1SBarry Smith .  text - short string describing option
138253acd3b1SBarry Smith .  man - manual page for option
138353acd3b1SBarry Smith -  nmax - maximum number of values
138453acd3b1SBarry Smith 
138553acd3b1SBarry Smith    Output Parameter:
138653acd3b1SBarry Smith +  value - location to copy values
138753acd3b1SBarry Smith .  nmax - actual number of values found
138853acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
138953acd3b1SBarry Smith 
139053acd3b1SBarry Smith    Level: beginner
139153acd3b1SBarry Smith 
139253acd3b1SBarry Smith    Notes:
139353acd3b1SBarry Smith    The user should pass in an array of doubles
139453acd3b1SBarry Smith 
139553acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
139653acd3b1SBarry Smith 
139753acd3b1SBarry Smith    Concepts: options database^array of strings
139853acd3b1SBarry Smith 
1399c5929fdfSBarry Smith .seealso: PetscOptionsGetInt(NULL,), PetscOptionsGetReal(),
1400acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
140153acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
140253acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1403acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1404a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
140553acd3b1SBarry Smith @*/
14064416b707SBarry Smith PetscErrorCode PetscOptionsRealArray_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscReal value[],PetscInt *n,PetscBool  *set)
140753acd3b1SBarry Smith {
140853acd3b1SBarry Smith   PetscErrorCode  ierr;
140953acd3b1SBarry Smith   PetscInt        i;
14104416b707SBarry Smith   PetscOptionItem amsopt;
141153acd3b1SBarry Smith 
141253acd3b1SBarry Smith   PetscFunctionBegin;
1413e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
1414e26ddf31SBarry Smith     PetscReal *vals;
1415e26ddf31SBarry Smith 
14164416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_REAL_ARRAY,&amsopt);CHKERRQ(ierr);
1417e55864a3SBarry Smith     ierr = PetscMalloc((*n)*sizeof(PetscReal),&amsopt->data);CHKERRQ(ierr);
1418e26ddf31SBarry Smith     vals = (PetscReal*)amsopt->data;
1419e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1420e26ddf31SBarry Smith     amsopt->arraylength = *n;
1421e26ddf31SBarry Smith   }
1422c5929fdfSBarry Smith   ierr = PetscOptionsGetRealArray(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,n,set);CHKERRQ(ierr);
1423e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1424a519f713SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%g",PetscOptionsObject->prefix?PetscOptionsObject->prefix:"",opt+1,(double)value[0]);CHKERRQ(ierr);
142553acd3b1SBarry Smith     for (i=1; i<*n; i++) {
1426e55864a3SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%g",(double)value[i]);CHKERRQ(ierr);
142753acd3b1SBarry Smith     }
1428e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
142953acd3b1SBarry Smith   }
143053acd3b1SBarry Smith   PetscFunctionReturn(0);
143153acd3b1SBarry Smith }
143253acd3b1SBarry Smith 
1433050cccc3SHong Zhang #undef __FUNCT__
1434050cccc3SHong Zhang #define __FUNCT__ "PetscOptionsScalarArray_Private"
1435050cccc3SHong Zhang /*@C
1436050cccc3SHong Zhang    PetscOptionsScalarArray - Gets an array of Scalar values for a particular
1437050cccc3SHong Zhang    option in the database. The values must be separated with commas with
1438050cccc3SHong Zhang    no intervening spaces.
1439050cccc3SHong Zhang 
1440050cccc3SHong Zhang    Logically Collective on the communicator passed in PetscOptionsBegin()
1441050cccc3SHong Zhang 
1442050cccc3SHong Zhang    Input Parameters:
1443050cccc3SHong Zhang +  opt - the option one is seeking
1444050cccc3SHong Zhang .  text - short string describing option
1445050cccc3SHong Zhang .  man - manual page for option
1446050cccc3SHong Zhang -  nmax - maximum number of values
1447050cccc3SHong Zhang 
1448050cccc3SHong Zhang    Output Parameter:
1449050cccc3SHong Zhang +  value - location to copy values
1450050cccc3SHong Zhang .  nmax - actual number of values found
1451050cccc3SHong Zhang -  set - PETSC_TRUE if found, else PETSC_FALSE
1452050cccc3SHong Zhang 
1453050cccc3SHong Zhang    Level: beginner
1454050cccc3SHong Zhang 
1455050cccc3SHong Zhang    Notes:
1456050cccc3SHong Zhang    The user should pass in an array of doubles
1457050cccc3SHong Zhang 
1458050cccc3SHong Zhang    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
1459050cccc3SHong Zhang 
1460050cccc3SHong Zhang    Concepts: options database^array of strings
1461050cccc3SHong Zhang 
1462c5929fdfSBarry Smith .seealso: PetscOptionsGetInt(NULL,), PetscOptionsGetReal(),
1463050cccc3SHong Zhang           PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
1464050cccc3SHong Zhang           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1465050cccc3SHong Zhang           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1466050cccc3SHong Zhang           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1467050cccc3SHong Zhang           PetscOptionsFList(), PetscOptionsEList()
1468050cccc3SHong Zhang @*/
14694416b707SBarry Smith PetscErrorCode PetscOptionsScalarArray_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscScalar value[],PetscInt *n,PetscBool  *set)
1470050cccc3SHong Zhang {
1471050cccc3SHong Zhang   PetscErrorCode  ierr;
1472050cccc3SHong Zhang   PetscInt        i;
14734416b707SBarry Smith   PetscOptionItem amsopt;
1474050cccc3SHong Zhang 
1475050cccc3SHong Zhang   PetscFunctionBegin;
1476050cccc3SHong Zhang   if (!PetscOptionsObject->count) {
1477050cccc3SHong Zhang     PetscScalar *vals;
1478050cccc3SHong Zhang 
14794416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_SCALAR_ARRAY,&amsopt);CHKERRQ(ierr);
1480050cccc3SHong Zhang     ierr = PetscMalloc((*n)*sizeof(PetscScalar),&amsopt->data);CHKERRQ(ierr);
1481050cccc3SHong Zhang     vals = (PetscScalar*)amsopt->data;
1482050cccc3SHong Zhang     for (i=0; i<*n; i++) vals[i] = value[i];
1483050cccc3SHong Zhang     amsopt->arraylength = *n;
1484050cccc3SHong Zhang   }
1485c5929fdfSBarry Smith   ierr = PetscOptionsGetScalarArray(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,n,set);CHKERRQ(ierr);
1486050cccc3SHong Zhang   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
148795f3a755SJose E. Roman     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%g+%gi",PetscOptionsObject->prefix?PetscOptionsObject->prefix:"",opt+1,(double)PetscRealPart(value[0]),(double)PetscImaginaryPart(value[0]));CHKERRQ(ierr);
1488050cccc3SHong Zhang     for (i=1; i<*n; i++) {
148995f3a755SJose E. Roman       ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%g+%gi",(double)PetscRealPart(value[i]),(double)PetscImaginaryPart(value[i]));CHKERRQ(ierr);
1490050cccc3SHong Zhang     }
1491050cccc3SHong Zhang     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
1492050cccc3SHong Zhang   }
1493050cccc3SHong Zhang   PetscFunctionReturn(0);
1494050cccc3SHong Zhang }
149553acd3b1SBarry Smith 
149653acd3b1SBarry Smith #undef __FUNCT__
1497e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsIntArray_Private"
149853acd3b1SBarry Smith /*@C
149953acd3b1SBarry Smith    PetscOptionsIntArray - Gets an array of integers for a particular
1500b32a342fSShri Abhyankar    option in the database.
150153acd3b1SBarry Smith 
15023f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
150353acd3b1SBarry Smith 
150453acd3b1SBarry Smith    Input Parameters:
150553acd3b1SBarry Smith +  opt - the option one is seeking
150653acd3b1SBarry Smith .  text - short string describing option
150753acd3b1SBarry Smith .  man - manual page for option
1508f8a50e2bSBarry Smith -  n - maximum number of values
150953acd3b1SBarry Smith 
151053acd3b1SBarry Smith    Output Parameter:
151153acd3b1SBarry Smith +  value - location to copy values
1512f8a50e2bSBarry Smith .  n - actual number of values found
151353acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
151453acd3b1SBarry Smith 
151553acd3b1SBarry Smith    Level: beginner
151653acd3b1SBarry Smith 
151753acd3b1SBarry Smith    Notes:
1518b32a342fSShri Abhyankar    The array can be passed as
1519bebe2cf6SSatish Balay    a comma separated list:                                 0,1,2,3,4,5,6,7
15200fd488f5SShri Abhyankar    a range (start-end+1):                                  0-8
15210fd488f5SShri Abhyankar    a range with given increment (start-end+1:inc):         0-7:2
1522bebe2cf6SSatish Balay    a combination of values and ranges separated by commas: 0,1-8,8-15:2
1523b32a342fSShri Abhyankar 
1524b32a342fSShri Abhyankar    There must be no intervening spaces between the values.
152553acd3b1SBarry Smith 
152653acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
152753acd3b1SBarry Smith 
1528b32a342fSShri Abhyankar    Concepts: options database^array of ints
152953acd3b1SBarry Smith 
1530c5929fdfSBarry Smith .seealso: PetscOptionsGetInt(NULL,), PetscOptionsGetReal(),
1531acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
153253acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
153353acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1534acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1535a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsRealArray()
153653acd3b1SBarry Smith @*/
15374416b707SBarry Smith PetscErrorCode  PetscOptionsIntArray_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscInt value[],PetscInt *n,PetscBool  *set)
153853acd3b1SBarry Smith {
153953acd3b1SBarry Smith   PetscErrorCode ierr;
154053acd3b1SBarry Smith   PetscInt        i;
15414416b707SBarry Smith   PetscOptionItem amsopt;
154253acd3b1SBarry Smith 
154353acd3b1SBarry Smith   PetscFunctionBegin;
1544e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
1545e26ddf31SBarry Smith     PetscInt *vals;
1546e26ddf31SBarry Smith 
15474416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_INT_ARRAY,&amsopt);CHKERRQ(ierr);
1548854ce69bSBarry Smith     ierr = PetscMalloc1(*n,(PetscInt**)&amsopt->data);CHKERRQ(ierr);
1549e26ddf31SBarry Smith     vals = (PetscInt*)amsopt->data;
1550e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1551e26ddf31SBarry Smith     amsopt->arraylength = *n;
1552e26ddf31SBarry Smith   }
1553c5929fdfSBarry Smith   ierr = PetscOptionsGetIntArray(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,n,set);CHKERRQ(ierr);
1554e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1555e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%d",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,value[0]);CHKERRQ(ierr);
155653acd3b1SBarry Smith     for (i=1; i<*n; i++) {
1557e55864a3SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%d",value[i]);CHKERRQ(ierr);
155853acd3b1SBarry Smith     }
1559e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
156053acd3b1SBarry Smith   }
156153acd3b1SBarry Smith   PetscFunctionReturn(0);
156253acd3b1SBarry Smith }
156353acd3b1SBarry Smith 
156453acd3b1SBarry Smith #undef __FUNCT__
1565e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsStringArray_Private"
156653acd3b1SBarry Smith /*@C
156753acd3b1SBarry Smith    PetscOptionsStringArray - Gets an array of string values for a particular
156853acd3b1SBarry Smith    option in the database. The values must be separated with commas with
156953acd3b1SBarry Smith    no intervening spaces.
157053acd3b1SBarry Smith 
15713f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
157253acd3b1SBarry Smith 
157353acd3b1SBarry Smith    Input Parameters:
157453acd3b1SBarry Smith +  opt - the option one is seeking
157553acd3b1SBarry Smith .  text - short string describing option
157653acd3b1SBarry Smith .  man - manual page for option
157753acd3b1SBarry Smith -  nmax - maximum number of strings
157853acd3b1SBarry Smith 
157953acd3b1SBarry Smith    Output Parameter:
158053acd3b1SBarry Smith +  value - location to copy strings
158153acd3b1SBarry Smith .  nmax - actual number of strings found
158253acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
158353acd3b1SBarry Smith 
158453acd3b1SBarry Smith    Level: beginner
158553acd3b1SBarry Smith 
158653acd3b1SBarry Smith    Notes:
158753acd3b1SBarry Smith    The user should pass in an array of pointers to char, to hold all the
158853acd3b1SBarry Smith    strings returned by this function.
158953acd3b1SBarry Smith 
159053acd3b1SBarry Smith    The user is responsible for deallocating the strings that are
159153acd3b1SBarry Smith    returned. The Fortran interface for this routine is not supported.
159253acd3b1SBarry Smith 
159353acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
159453acd3b1SBarry Smith 
159553acd3b1SBarry Smith    Concepts: options database^array of strings
159653acd3b1SBarry Smith 
1597c5929fdfSBarry Smith .seealso: PetscOptionsGetInt(NULL,), PetscOptionsGetReal(),
1598acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
159953acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
160053acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1601acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1602a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
160353acd3b1SBarry Smith @*/
16044416b707SBarry Smith PetscErrorCode  PetscOptionsStringArray_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],char *value[],PetscInt *nmax,PetscBool  *set)
160553acd3b1SBarry Smith {
160653acd3b1SBarry Smith   PetscErrorCode  ierr;
16074416b707SBarry Smith   PetscOptionItem amsopt;
160853acd3b1SBarry Smith 
160953acd3b1SBarry Smith   PetscFunctionBegin;
1610e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
16114416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_STRING_ARRAY,&amsopt);CHKERRQ(ierr);
1612854ce69bSBarry Smith     ierr = PetscMalloc1(*nmax,(char**)&amsopt->data);CHKERRQ(ierr);
1613a297a907SKarl Rupp 
16141ae3d29cSBarry Smith     amsopt->arraylength = *nmax;
16151ae3d29cSBarry Smith   }
1616c5929fdfSBarry Smith   ierr = PetscOptionsGetStringArray(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,nmax,set);CHKERRQ(ierr);
1617e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1618e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <string1,string2,...>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
161953acd3b1SBarry Smith   }
162053acd3b1SBarry Smith   PetscFunctionReturn(0);
162153acd3b1SBarry Smith }
162253acd3b1SBarry Smith 
1623e2446a98SMatthew Knepley #undef __FUNCT__
1624e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsBoolArray_Private"
1625e2446a98SMatthew Knepley /*@C
1626acfcf0e5SJed Brown    PetscOptionsBoolArray - Gets an array of logical values (true or false) for a particular
1627e2446a98SMatthew Knepley    option in the database. The values must be separated with commas with
1628e2446a98SMatthew Knepley    no intervening spaces.
1629e2446a98SMatthew Knepley 
16303f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
1631e2446a98SMatthew Knepley 
1632e2446a98SMatthew Knepley    Input Parameters:
1633e2446a98SMatthew Knepley +  opt - the option one is seeking
1634e2446a98SMatthew Knepley .  text - short string describing option
1635e2446a98SMatthew Knepley .  man - manual page for option
1636e2446a98SMatthew Knepley -  nmax - maximum number of values
1637e2446a98SMatthew Knepley 
1638e2446a98SMatthew Knepley    Output Parameter:
1639e2446a98SMatthew Knepley +  value - location to copy values
1640e2446a98SMatthew Knepley .  nmax - actual number of values found
1641e2446a98SMatthew Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
1642e2446a98SMatthew Knepley 
1643e2446a98SMatthew Knepley    Level: beginner
1644e2446a98SMatthew Knepley 
1645e2446a98SMatthew Knepley    Notes:
1646e2446a98SMatthew Knepley    The user should pass in an array of doubles
1647e2446a98SMatthew Knepley 
1648e2446a98SMatthew Knepley    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
1649e2446a98SMatthew Knepley 
1650e2446a98SMatthew Knepley    Concepts: options database^array of strings
1651e2446a98SMatthew Knepley 
1652c5929fdfSBarry Smith .seealso: PetscOptionsGetInt(NULL,), PetscOptionsGetReal(),
1653acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
1654e2446a98SMatthew Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1655e2446a98SMatthew Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1656acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1657a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
1658e2446a98SMatthew Knepley @*/
16594416b707SBarry Smith PetscErrorCode  PetscOptionsBoolArray_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool value[],PetscInt *n,PetscBool *set)
1660e2446a98SMatthew Knepley {
1661e2446a98SMatthew Knepley   PetscErrorCode   ierr;
1662e2446a98SMatthew Knepley   PetscInt         i;
16634416b707SBarry Smith   PetscOptionItem  amsopt;
1664e2446a98SMatthew Knepley 
1665e2446a98SMatthew Knepley   PetscFunctionBegin;
1666e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
1667ace3abfcSBarry Smith     PetscBool *vals;
16681ae3d29cSBarry Smith 
16694416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL_ARRAY,&amsopt);CHKERRQ(ierr);
16701a1499c8SBarry Smith     ierr = PetscMalloc1(*n,(PetscBool**)&amsopt->data);CHKERRQ(ierr);
1671ace3abfcSBarry Smith     vals = (PetscBool*)amsopt->data;
16721ae3d29cSBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
16731ae3d29cSBarry Smith     amsopt->arraylength = *n;
16741ae3d29cSBarry Smith   }
1675c5929fdfSBarry Smith   ierr = PetscOptionsGetBoolArray(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,n,set);CHKERRQ(ierr);
1676e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1677e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%d",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,value[0]);CHKERRQ(ierr);
1678e2446a98SMatthew Knepley     for (i=1; i<*n; i++) {
1679e55864a3SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%d",value[i]);CHKERRQ(ierr);
1680e2446a98SMatthew Knepley     }
1681e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
1682e2446a98SMatthew Knepley   }
1683e2446a98SMatthew Knepley   PetscFunctionReturn(0);
1684e2446a98SMatthew Knepley }
1685e2446a98SMatthew Knepley 
16868cc676e6SMatthew G Knepley #undef __FUNCT__
1687e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsViewer_Private"
16888cc676e6SMatthew G Knepley /*@C
1689d1da0b69SBarry Smith    PetscOptionsViewer - Gets a viewer appropriate for the type indicated by the user
16908cc676e6SMatthew G Knepley 
16918cc676e6SMatthew G Knepley    Logically Collective on the communicator passed in PetscOptionsBegin()
16928cc676e6SMatthew G Knepley 
16938cc676e6SMatthew G Knepley    Input Parameters:
16948cc676e6SMatthew G Knepley +  opt - option name
16958cc676e6SMatthew G Knepley .  text - short string that describes the option
16968cc676e6SMatthew G Knepley -  man - manual page with additional information on option
16978cc676e6SMatthew G Knepley 
16988cc676e6SMatthew G Knepley    Output Parameter:
16998cc676e6SMatthew G Knepley +  viewer - the viewer
17008cc676e6SMatthew G Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
17018cc676e6SMatthew G Knepley 
17028cc676e6SMatthew G Knepley    Level: beginner
17038cc676e6SMatthew G Knepley 
17048cc676e6SMatthew G Knepley    Concepts: options database^has int
17058cc676e6SMatthew G Knepley 
17068cc676e6SMatthew G Knepley    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
17078cc676e6SMatthew G Knepley 
17085a7113b9SPatrick Sanan    See PetscOptionsGetViewer() for the format of the supplied viewer and its options
17098cc676e6SMatthew G Knepley 
1710c5929fdfSBarry Smith .seealso: PetscOptionsGetViewer(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(NULL,),
17118cc676e6SMatthew G Knepley           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
17128cc676e6SMatthew G Knepley           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
17138cc676e6SMatthew G Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
17148cc676e6SMatthew G Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
17158cc676e6SMatthew G Knepley           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1716a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
17178cc676e6SMatthew G Knepley @*/
17184416b707SBarry Smith PetscErrorCode  PetscOptionsViewer_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscViewer *viewer,PetscViewerFormat *format,PetscBool  *set)
17198cc676e6SMatthew G Knepley {
17208cc676e6SMatthew G Knepley   PetscErrorCode  ierr;
17214416b707SBarry Smith   PetscOptionItem amsopt;
17228cc676e6SMatthew G Knepley 
17238cc676e6SMatthew G Knepley   PetscFunctionBegin;
17241a1499c8SBarry Smith   if (!PetscOptionsObject->count) {
17254416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
172664facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
17275b02f95dSBarry Smith     ierr = PetscStrdup("",(char**)&amsopt->data);CHKERRQ(ierr);
17288cc676e6SMatthew G Knepley   }
1729e55864a3SBarry Smith   ierr = PetscOptionsGetViewer(PetscOptionsObject->comm,PetscOptionsObject->prefix,opt,viewer,format,set);CHKERRQ(ierr);
1730e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1731e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%s>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,"",text,ManSection(man));CHKERRQ(ierr);
17328cc676e6SMatthew G Knepley   }
17338cc676e6SMatthew G Knepley   PetscFunctionReturn(0);
17348cc676e6SMatthew G Knepley }
17358cc676e6SMatthew G Knepley 
173653acd3b1SBarry Smith 
173753acd3b1SBarry Smith #undef __FUNCT__
173853acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsHead"
173953acd3b1SBarry Smith /*@C
1740b52f573bSBarry Smith      PetscOptionsHead - Puts a heading before listing any more published options. Used, for example,
174153acd3b1SBarry Smith             in KSPSetFromOptions_GMRES().
174253acd3b1SBarry Smith 
17433f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
174453acd3b1SBarry Smith 
174553acd3b1SBarry Smith    Input Parameter:
174653acd3b1SBarry Smith .   head - the heading text
174753acd3b1SBarry Smith 
174853acd3b1SBarry Smith 
174953acd3b1SBarry Smith    Level: intermediate
175053acd3b1SBarry Smith 
175153acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
175253acd3b1SBarry Smith 
1753b52f573bSBarry Smith           Can be followed by a call to PetscOptionsTail() in the same function.
175453acd3b1SBarry Smith 
175553acd3b1SBarry Smith    Concepts: options database^subheading
175653acd3b1SBarry Smith 
1757c5929fdfSBarry Smith .seealso: PetscOptionsGetInt(NULL,), PetscOptionsGetReal(),
1758acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
175953acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
176053acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1761acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1762a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
176353acd3b1SBarry Smith @*/
17644416b707SBarry Smith PetscErrorCode  PetscOptionsHead(PetscOptionItems *PetscOptionsObject,const char head[])
176553acd3b1SBarry Smith {
176653acd3b1SBarry Smith   PetscErrorCode ierr;
176753acd3b1SBarry Smith 
176853acd3b1SBarry Smith   PetscFunctionBegin;
1769e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1770e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  %s\n",head);CHKERRQ(ierr);
177153acd3b1SBarry Smith   }
177253acd3b1SBarry Smith   PetscFunctionReturn(0);
177353acd3b1SBarry Smith }
177453acd3b1SBarry Smith 
177553acd3b1SBarry Smith 
177653acd3b1SBarry Smith 
177753acd3b1SBarry Smith 
177853acd3b1SBarry Smith 
177953acd3b1SBarry Smith 
1780