xref: /petsc/src/sys/objects/aoptions.c (revision 570b7f6d312edbd2962258513e9546951e9f0a9d)
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;
320eb63584SBarry Smith   if (!PetscOptionsObject->alreadyprinted) {
339de0f6ecSBarry Smith     if (!PetscOptionsHelpPrintedSingleton) {
349de0f6ecSBarry Smith       ierr = PetscOptionsHelpPrintedCreate(&PetscOptionsHelpPrintedSingleton);CHKERRQ(ierr);
359de0f6ecSBarry Smith     }
369de0f6ecSBarry Smith     ierr = PetscOptionsHelpPrintedCheck(PetscOptionsHelpPrintedSingleton,prefix,title,&PetscOptionsObject->alreadyprinted);CHKERRQ(ierr);
370eb63584SBarry Smith   }
38e55864a3SBarry Smith   PetscOptionsObject->next          = 0;
39e55864a3SBarry Smith   PetscOptionsObject->comm          = comm;
40e55864a3SBarry Smith   PetscOptionsObject->changedmethod = PETSC_FALSE;
41a297a907SKarl Rupp 
42e55864a3SBarry Smith   ierr = PetscStrallocpy(prefix,&PetscOptionsObject->prefix);CHKERRQ(ierr);
43e55864a3SBarry Smith   ierr = PetscStrallocpy(title,&PetscOptionsObject->title);CHKERRQ(ierr);
4453acd3b1SBarry Smith 
45c5929fdfSBarry Smith   ierr = PetscOptionsHasName(PetscOptionsObject->options,NULL,"-help",&PetscOptionsObject->printhelp);CHKERRQ(ierr);
46e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1) {
47e55864a3SBarry Smith     if (!PetscOptionsObject->alreadyprinted) {
4853acd3b1SBarry Smith       ierr = (*PetscHelpPrintf)(comm,"%s -------------------------------------------------\n",title);CHKERRQ(ierr);
4953acd3b1SBarry Smith     }
5061b37b28SSatish Balay   }
5153acd3b1SBarry Smith   PetscFunctionReturn(0);
5253acd3b1SBarry Smith }
5353acd3b1SBarry Smith 
543194b578SJed Brown #undef __FUNCT__
553194b578SJed Brown #define __FUNCT__ "PetscObjectOptionsBegin_Private"
563194b578SJed Brown /*
573194b578SJed Brown     Handles setting up the data structure in a call to PetscObjectOptionsBegin()
583194b578SJed Brown */
594416b707SBarry Smith PetscErrorCode PetscObjectOptionsBegin_Private(PetscOptionItems *PetscOptionsObject,PetscObject obj)
603194b578SJed Brown {
613194b578SJed Brown   PetscErrorCode ierr;
623194b578SJed Brown   char           title[256];
633194b578SJed Brown   PetscBool      flg;
643194b578SJed Brown 
653194b578SJed Brown   PetscFunctionBegin;
663194b578SJed Brown   PetscValidHeader(obj,1);
67e55864a3SBarry Smith   PetscOptionsObject->object         = obj;
68e55864a3SBarry Smith   PetscOptionsObject->alreadyprinted = obj->optionsprinted;
69a297a907SKarl Rupp 
703194b578SJed Brown   ierr = PetscStrcmp(obj->description,obj->class_name,&flg);CHKERRQ(ierr);
713194b578SJed Brown   if (flg) {
728caf3d72SBarry Smith     ierr = PetscSNPrintf(title,sizeof(title),"%s options",obj->class_name);CHKERRQ(ierr);
733194b578SJed Brown   } else {
748caf3d72SBarry Smith     ierr = PetscSNPrintf(title,sizeof(title),"%s (%s) options",obj->description,obj->class_name);CHKERRQ(ierr);
753194b578SJed Brown   }
76e55864a3SBarry Smith   ierr = PetscOptionsBegin_Private(PetscOptionsObject,obj->comm,obj->prefix,title,obj->mansec);CHKERRQ(ierr);
773194b578SJed Brown   PetscFunctionReturn(0);
783194b578SJed Brown }
793194b578SJed Brown 
8053acd3b1SBarry Smith /*
8153acd3b1SBarry Smith      Handles adding another option to the list of options within this particular PetscOptionsBegin() PetscOptionsEnd()
8253acd3b1SBarry Smith */
8353acd3b1SBarry Smith #undef __FUNCT__
848a19d54fSBarry Smith #define __FUNCT__ "PetscOptionItemCreate_Private"
854416b707SBarry Smith static int PetscOptionItemCreate_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscOptionType t,PetscOptionItem *amsopt)
8653acd3b1SBarry Smith {
8753acd3b1SBarry Smith   int             ierr;
884416b707SBarry Smith   PetscOptionItem next;
893be6e4c3SJed Brown   PetscBool       valid;
9053acd3b1SBarry Smith 
9153acd3b1SBarry Smith   PetscFunctionBegin;
923be6e4c3SJed Brown   ierr = PetscOptionsValidKey(opt,&valid);CHKERRQ(ierr);
933be6e4c3SJed Brown   if (!valid) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"The option '%s' is not a valid key",opt);
943be6e4c3SJed Brown 
95b00a9115SJed Brown   ierr            = PetscNew(amsopt);CHKERRQ(ierr);
9653acd3b1SBarry Smith   (*amsopt)->next = 0;
9753acd3b1SBarry Smith   (*amsopt)->set  = PETSC_FALSE;
986356e834SBarry Smith   (*amsopt)->type = t;
9953acd3b1SBarry Smith   (*amsopt)->data = 0;
10061b37b28SSatish Balay 
10153acd3b1SBarry Smith   ierr = PetscStrallocpy(text,&(*amsopt)->text);CHKERRQ(ierr);
10253acd3b1SBarry Smith   ierr = PetscStrallocpy(opt,&(*amsopt)->option);CHKERRQ(ierr);
1036356e834SBarry Smith   ierr = PetscStrallocpy(man,&(*amsopt)->man);CHKERRQ(ierr);
10453acd3b1SBarry Smith 
105e55864a3SBarry Smith   if (!PetscOptionsObject->next) PetscOptionsObject->next = *amsopt;
106a297a907SKarl Rupp   else {
107e55864a3SBarry Smith     next = PetscOptionsObject->next;
10853acd3b1SBarry Smith     while (next->next) next = next->next;
10953acd3b1SBarry Smith     next->next = *amsopt;
11053acd3b1SBarry Smith   }
11153acd3b1SBarry Smith   PetscFunctionReturn(0);
11253acd3b1SBarry Smith }
11353acd3b1SBarry Smith 
11453acd3b1SBarry Smith #undef __FUNCT__
115aee2cecaSBarry Smith #define __FUNCT__ "PetscScanString"
116aee2cecaSBarry Smith /*
1173fc1eb6aSBarry Smith     PetscScanString -  Gets user input via stdin from process and broadcasts to all processes
1183fc1eb6aSBarry Smith 
1193fc1eb6aSBarry Smith     Collective on MPI_Comm
1203fc1eb6aSBarry Smith 
1213fc1eb6aSBarry Smith    Input Parameters:
1223fc1eb6aSBarry Smith +     commm - communicator for the broadcast, must be PETSC_COMM_WORLD
1233fc1eb6aSBarry Smith .     n - length of the string, must be the same on all processes
1243fc1eb6aSBarry Smith -     str - location to store input
125aee2cecaSBarry Smith 
126aee2cecaSBarry Smith     Bugs:
127aee2cecaSBarry Smith .   Assumes process 0 of the given communicator has access to stdin
128aee2cecaSBarry Smith 
129aee2cecaSBarry Smith */
1303fc1eb6aSBarry Smith static PetscErrorCode PetscScanString(MPI_Comm comm,size_t n,char str[])
131aee2cecaSBarry Smith {
132330cf3c9SBarry Smith   size_t         i;
133aee2cecaSBarry Smith   char           c;
1343fc1eb6aSBarry Smith   PetscMPIInt    rank,nm;
135aee2cecaSBarry Smith   PetscErrorCode ierr;
136aee2cecaSBarry Smith 
137aee2cecaSBarry Smith   PetscFunctionBegin;
138aee2cecaSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
139aee2cecaSBarry Smith   if (!rank) {
140aee2cecaSBarry Smith     c = (char) getchar();
141aee2cecaSBarry Smith     i = 0;
142aee2cecaSBarry Smith     while (c != '\n' && i < n-1) {
143aee2cecaSBarry Smith       str[i++] = c;
144aee2cecaSBarry Smith       c = (char)getchar();
145aee2cecaSBarry Smith     }
146aee2cecaSBarry Smith     str[i] = 0;
147aee2cecaSBarry Smith   }
1484dc2109aSBarry Smith   ierr = PetscMPIIntCast(n,&nm);CHKERRQ(ierr);
1493fc1eb6aSBarry Smith   ierr = MPI_Bcast(str,nm,MPI_CHAR,0,comm);CHKERRQ(ierr);
150aee2cecaSBarry Smith   PetscFunctionReturn(0);
151aee2cecaSBarry Smith }
152aee2cecaSBarry Smith 
153ead66b60SBarry Smith #undef __FUNCT__
154ead66b60SBarry Smith #define __FUNCT__ "PetscStrdup"
1555b02f95dSBarry Smith /*
1565b02f95dSBarry Smith     This is needed because certain strings may be freed by SAWs, hence we cannot use PetscStrallocpy()
1575b02f95dSBarry Smith */
1585b02f95dSBarry Smith static PetscErrorCode  PetscStrdup(const char s[],char *t[])
1595b02f95dSBarry Smith {
1605b02f95dSBarry Smith   PetscErrorCode ierr;
1615b02f95dSBarry Smith   size_t         len;
1625b02f95dSBarry Smith   char           *tmp = 0;
1635b02f95dSBarry Smith 
1645b02f95dSBarry Smith   PetscFunctionBegin;
1655b02f95dSBarry Smith   if (s) {
1665b02f95dSBarry Smith     ierr = PetscStrlen(s,&len);CHKERRQ(ierr);
167f416af30SBarry Smith     tmp = (char*) malloc((len+1)*sizeof(char));
1685b02f95dSBarry Smith     if (!tmp) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"No memory to duplicate string");
1695b02f95dSBarry Smith     ierr = PetscStrcpy(tmp,s);CHKERRQ(ierr);
1705b02f95dSBarry Smith   }
1715b02f95dSBarry Smith   *t = tmp;
1725b02f95dSBarry Smith   PetscFunctionReturn(0);
1735b02f95dSBarry Smith }
1745b02f95dSBarry Smith 
1755b02f95dSBarry Smith 
176aee2cecaSBarry Smith #undef __FUNCT__
177aee2cecaSBarry Smith #define __FUNCT__ "PetscOptionsGetFromTextInput"
178aee2cecaSBarry Smith /*
1793cc1e11dSBarry Smith     PetscOptionsGetFromTextInput - Presents all the PETSc Options processed by the program so the user may change them at runtime
180aee2cecaSBarry Smith 
181aee2cecaSBarry Smith     Notes: this isn't really practical, it is just to demonstrate the principle
182aee2cecaSBarry Smith 
1837781c08eSBarry Smith     A carriage return indicates no change from the default; but this like -ksp_monitor <stdout>  the default is actually not stdout the default
1847781c08eSBarry Smith     is to do nothing so to get it to use stdout you need to type stdout. This is kind of bug?
1857781c08eSBarry Smith 
186aee2cecaSBarry Smith     Bugs:
1877781c08eSBarry Smith +    All processes must traverse through the exact same set of option queries due to the call to PetscScanString()
1883cc1e11dSBarry Smith .    Internal strings have arbitrary length and string copies are not checked that they fit into string space
189aee2cecaSBarry Smith -    Only works for PetscInt == int, PetscReal == double etc
190aee2cecaSBarry Smith 
1913cc1e11dSBarry Smith     Developer Notes: Normally the GUI that presents the options the user and retrieves the values would be running in a different
1923cc1e11dSBarry Smith      address space and communicating with the PETSc program
1933cc1e11dSBarry Smith 
194aee2cecaSBarry Smith */
1954416b707SBarry Smith PetscErrorCode PetscOptionsGetFromTextInput(PetscOptionItems *PetscOptionsObject)
1966356e834SBarry Smith {
1976356e834SBarry Smith   PetscErrorCode  ierr;
1984416b707SBarry Smith   PetscOptionItem next = PetscOptionsObject->next;
1996356e834SBarry Smith   char            str[512];
2007781c08eSBarry Smith   PetscBool       bid;
201a4404d99SBarry Smith   PetscReal       ir,*valr;
202330cf3c9SBarry Smith   PetscInt        *vald;
203330cf3c9SBarry Smith   size_t          i;
2046356e834SBarry Smith 
2052a409bb0SBarry Smith   PetscFunctionBegin;
206e55864a3SBarry Smith   ierr = (*PetscPrintf)(PETSC_COMM_WORLD,"%s -------------------------------------------------\n",PetscOptionsObject->title);CHKERRQ(ierr);
2076356e834SBarry Smith   while (next) {
2086356e834SBarry Smith     switch (next->type) {
2096356e834SBarry Smith     case OPTION_HEAD:
2106356e834SBarry Smith       break;
211e26ddf31SBarry Smith     case OPTION_INT_ARRAY:
212e55864a3SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",next->option+1);CHKERRQ(ierr);
213e26ddf31SBarry Smith       vald = (PetscInt*) next->data;
214e26ddf31SBarry Smith       for (i=0; i<next->arraylength; i++) {
215e26ddf31SBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"%d",vald[i]);CHKERRQ(ierr);
216e26ddf31SBarry Smith         if (i < next->arraylength-1) {
217e26ddf31SBarry Smith           ierr = PetscPrintf(PETSC_COMM_WORLD,",");CHKERRQ(ierr);
218e26ddf31SBarry Smith         }
219e26ddf31SBarry Smith       }
220e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,">: %s (%s) ",next->text,next->man);CHKERRQ(ierr);
221e26ddf31SBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
222e26ddf31SBarry Smith       if (str[0]) {
223e26ddf31SBarry Smith         PetscToken token;
224e26ddf31SBarry Smith         PetscInt   n=0,nmax = next->arraylength,*dvalue = (PetscInt*)next->data,start,end;
225e26ddf31SBarry Smith         size_t     len;
226e26ddf31SBarry Smith         char       *value;
227ace3abfcSBarry Smith         PetscBool  foundrange;
228e26ddf31SBarry Smith 
229e26ddf31SBarry Smith         next->set = PETSC_TRUE;
230e26ddf31SBarry Smith         value     = str;
231e26ddf31SBarry Smith         ierr      = PetscTokenCreate(value,',',&token);CHKERRQ(ierr);
232e26ddf31SBarry Smith         ierr      = PetscTokenFind(token,&value);CHKERRQ(ierr);
233e26ddf31SBarry Smith         while (n < nmax) {
234e26ddf31SBarry Smith           if (!value) break;
235e26ddf31SBarry Smith 
236e26ddf31SBarry Smith           /* look for form  d-D where d and D are integers */
237e26ddf31SBarry Smith           foundrange = PETSC_FALSE;
238e26ddf31SBarry Smith           ierr       = PetscStrlen(value,&len);CHKERRQ(ierr);
239e26ddf31SBarry Smith           if (value[0] == '-') i=2;
240e26ddf31SBarry Smith           else i=1;
241330cf3c9SBarry Smith           for (;i<len; i++) {
242e26ddf31SBarry Smith             if (value[i] == '-') {
243e32f2f54SBarry Smith               if (i == len-1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Error in %D-th array entry %s\n",n,value);
244e26ddf31SBarry Smith               value[i] = 0;
245cfbddea1SSatish Balay               ierr     = PetscOptionsStringToInt(value,&start);CHKERRQ(ierr);
246cfbddea1SSatish Balay               ierr     = PetscOptionsStringToInt(value+i+1,&end);CHKERRQ(ierr);
247e32f2f54SBarry 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);
248e32f2f54SBarry 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);
249e26ddf31SBarry Smith               for (; start<end; start++) {
250e26ddf31SBarry Smith                 *dvalue = start; dvalue++;n++;
251e26ddf31SBarry Smith               }
252e26ddf31SBarry Smith               foundrange = PETSC_TRUE;
253e26ddf31SBarry Smith               break;
254e26ddf31SBarry Smith             }
255e26ddf31SBarry Smith           }
256e26ddf31SBarry Smith           if (!foundrange) {
257cfbddea1SSatish Balay             ierr = PetscOptionsStringToInt(value,dvalue);CHKERRQ(ierr);
258e26ddf31SBarry Smith             dvalue++;
259e26ddf31SBarry Smith             n++;
260e26ddf31SBarry Smith           }
261e26ddf31SBarry Smith           ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
262e26ddf31SBarry Smith         }
2638c74ee41SBarry Smith         ierr = PetscTokenDestroy(&token);CHKERRQ(ierr);
264e26ddf31SBarry Smith       }
265e26ddf31SBarry Smith       break;
266e26ddf31SBarry Smith     case OPTION_REAL_ARRAY:
267e55864a3SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",next->option+1);CHKERRQ(ierr);
268e26ddf31SBarry Smith       valr = (PetscReal*) next->data;
269e26ddf31SBarry Smith       for (i=0; i<next->arraylength; i++) {
270e26ddf31SBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"%g",valr[i]);CHKERRQ(ierr);
271e26ddf31SBarry Smith         if (i < next->arraylength-1) {
272e26ddf31SBarry Smith           ierr = PetscPrintf(PETSC_COMM_WORLD,",");CHKERRQ(ierr);
273e26ddf31SBarry Smith         }
274e26ddf31SBarry Smith       }
275e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,">: %s (%s) ",next->text,next->man);CHKERRQ(ierr);
276e26ddf31SBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
277e26ddf31SBarry Smith       if (str[0]) {
278e26ddf31SBarry Smith         PetscToken token;
279e26ddf31SBarry Smith         PetscInt   n = 0,nmax = next->arraylength;
280e26ddf31SBarry Smith         PetscReal  *dvalue = (PetscReal*)next->data;
281e26ddf31SBarry Smith         char       *value;
282e26ddf31SBarry Smith 
283e26ddf31SBarry Smith         next->set = PETSC_TRUE;
284e26ddf31SBarry Smith         value     = str;
285e26ddf31SBarry Smith         ierr      = PetscTokenCreate(value,',',&token);CHKERRQ(ierr);
286e26ddf31SBarry Smith         ierr      = PetscTokenFind(token,&value);CHKERRQ(ierr);
287e26ddf31SBarry Smith         while (n < nmax) {
288e26ddf31SBarry Smith           if (!value) break;
289cfbddea1SSatish Balay           ierr = PetscOptionsStringToReal(value,dvalue);CHKERRQ(ierr);
290e26ddf31SBarry Smith           dvalue++;
291e26ddf31SBarry Smith           n++;
292e26ddf31SBarry Smith           ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
293e26ddf31SBarry Smith         }
2948c74ee41SBarry Smith         ierr = PetscTokenDestroy(&token);CHKERRQ(ierr);
295e26ddf31SBarry Smith       }
296e26ddf31SBarry Smith       break;
2976356e834SBarry Smith     case OPTION_INT:
298e55864a3SBarry 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);
2993fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3003fc1eb6aSBarry Smith       if (str[0]) {
301d25d7f95SJed Brown #if defined(PETSC_SIZEOF_LONG_LONG)
302d25d7f95SJed Brown         long long lid;
303d25d7f95SJed Brown         sscanf(str,"%lld",&lid);
3041a1499c8SBarry 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);
305c272547aSJed Brown #else
306d25d7f95SJed Brown         long  lid;
307d25d7f95SJed Brown         sscanf(str,"%ld",&lid);
3081a1499c8SBarry 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);
309c272547aSJed Brown #endif
310a297a907SKarl Rupp 
311d25d7f95SJed Brown         next->set = PETSC_TRUE;
312d25d7f95SJed Brown         *((PetscInt*)next->data) = (PetscInt)lid;
313aee2cecaSBarry Smith       }
314aee2cecaSBarry Smith       break;
315aee2cecaSBarry Smith     case OPTION_REAL:
316e55864a3SBarry 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);
3173fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3183fc1eb6aSBarry Smith       if (str[0]) {
319ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
320a4404d99SBarry Smith         sscanf(str,"%e",&ir);
321*570b7f6dSBarry Smith #elif defined(PETSC_USE_REAL___FP16)
322*570b7f6dSBarry Smith         float irtemp;
323*570b7f6dSBarry Smith         sscanf(str,"%e",&irtemp);
324*570b7f6dSBarry Smith         ir = irtemp;
325ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL_DOUBLE)
326aee2cecaSBarry Smith         sscanf(str,"%le",&ir);
327ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL___FLOAT128)
328d9822059SBarry Smith         ir = strtoflt128(str,0);
329d9822059SBarry Smith #else
330513dbe71SLisandro Dalcin         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Unknown scalar type");
331a4404d99SBarry Smith #endif
332aee2cecaSBarry Smith         next->set                 = PETSC_TRUE;
333aee2cecaSBarry Smith         *((PetscReal*)next->data) = ir;
334aee2cecaSBarry Smith       }
335aee2cecaSBarry Smith       break;
3367781c08eSBarry Smith     case OPTION_BOOL:
33783355fc5SBarry 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);
3387781c08eSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3397781c08eSBarry Smith       if (str[0]) {
3407781c08eSBarry Smith         ierr = PetscOptionsStringToBool(str,&bid);CHKERRQ(ierr);
3417781c08eSBarry Smith         next->set = PETSC_TRUE;
3427781c08eSBarry Smith         *((PetscBool*)next->data) = bid;
3437781c08eSBarry Smith       }
3447781c08eSBarry Smith       break;
345aee2cecaSBarry Smith     case OPTION_STRING:
346e55864a3SBarry 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);
3473fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3483fc1eb6aSBarry Smith       if (str[0]) {
349aee2cecaSBarry Smith         next->set = PETSC_TRUE;
35064facd6cSBarry Smith         /* must use system malloc since SAWs may free this */
3515b02f95dSBarry Smith         ierr = PetscStrdup(str,(char**)&next->data);CHKERRQ(ierr);
3526356e834SBarry Smith       }
3536356e834SBarry Smith       break;
354a264d7a6SBarry Smith     case OPTION_FLIST:
355e55864a3SBarry Smith       ierr = PetscFunctionListPrintTypes(PETSC_COMM_WORLD,stdout,PetscOptionsObject->prefix,next->option,next->text,next->man,next->flist,(char*)next->data);CHKERRQ(ierr);
3563cc1e11dSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3573cc1e11dSBarry Smith       if (str[0]) {
358e55864a3SBarry Smith         PetscOptionsObject->changedmethod = PETSC_TRUE;
3593cc1e11dSBarry Smith         next->set = PETSC_TRUE;
36064facd6cSBarry Smith         /* must use system malloc since SAWs may free this */
3615b02f95dSBarry Smith         ierr = PetscStrdup(str,(char**)&next->data);CHKERRQ(ierr);
3623cc1e11dSBarry Smith       }
3633cc1e11dSBarry Smith       break;
364b432afdaSMatthew Knepley     default:
365b432afdaSMatthew Knepley       break;
3666356e834SBarry Smith     }
3676356e834SBarry Smith     next = next->next;
3686356e834SBarry Smith   }
3696356e834SBarry Smith   PetscFunctionReturn(0);
3706356e834SBarry Smith }
3716356e834SBarry Smith 
372e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
373e04113cfSBarry Smith #include <petscviewersaws.h>
374d5649816SBarry Smith 
375d5649816SBarry Smith static int count = 0;
376d5649816SBarry Smith 
377b3506946SBarry Smith #undef __FUNCT__
378e04113cfSBarry Smith #define __FUNCT__ "PetscOptionsSAWsDestroy"
379e04113cfSBarry Smith PetscErrorCode PetscOptionsSAWsDestroy(void)
380d5649816SBarry Smith {
3812657e9d9SBarry Smith   PetscFunctionBegin;
382d5649816SBarry Smith   PetscFunctionReturn(0);
383d5649816SBarry Smith }
384d5649816SBarry Smith 
3859c1e0ce8SBarry Smith static const char *OptionsHeader = "<head>\n"
38623a1ff79SBarry Smith                                    "<script type=\"text/javascript\" src=\"http://www.mcs.anl.gov/research/projects/saws/js/jquery-1.9.1.js\"></script>\n"
38723a1ff79SBarry Smith                                    "<script type=\"text/javascript\" src=\"http://www.mcs.anl.gov/research/projects/saws/js/SAWs.js\"></script>\n"
388d1fc0251SBarry Smith                                    "<script type=\"text/javascript\" src=\"js/PETSc.js\"></script>\n"
38964bbc9efSBarry Smith                                    "<script>\n"
39064bbc9efSBarry Smith                                       "jQuery(document).ready(function() {\n"
39164bbc9efSBarry Smith                                          "PETSc.getAndDisplayDirectory(null,\"#variablesInfo\")\n"
39264bbc9efSBarry Smith                                       "})\n"
39364bbc9efSBarry Smith                                   "</script>\n"
39464bbc9efSBarry Smith                                   "</head>\n";
3951423471aSBarry Smith 
3961423471aSBarry Smith /*  Determines the size and style of the scroll region where PETSc options selectable from users are displayed */
3971423471aSBarry Smith static const char *OptionsBodyBottom = "<div id=\"variablesInfo\" style=\"background-color:lightblue;height:auto;max-height:500px;overflow:scroll;\"></div>\n<br>\n</body>";
39864bbc9efSBarry Smith 
399d5649816SBarry Smith #undef __FUNCT__
4007781c08eSBarry Smith #define __FUNCT__ "PetscOptionsSAWsInput"
401b3506946SBarry Smith /*
4027781c08eSBarry Smith     PetscOptionsSAWsInput - Presents all the PETSc Options processed by the program so the user may change them at runtime using the SAWs
403b3506946SBarry Smith 
404b3506946SBarry Smith     Bugs:
405b3506946SBarry Smith +    All processes must traverse through the exact same set of option queries do to the call to PetscScanString()
406b3506946SBarry Smith .    Internal strings have arbitrary length and string copies are not checked that they fit into string space
407b3506946SBarry Smith -    Only works for PetscInt == int, PetscReal == double etc
408b3506946SBarry Smith 
409b3506946SBarry Smith 
410b3506946SBarry Smith */
4114416b707SBarry Smith PetscErrorCode PetscOptionsSAWsInput(PetscOptionItems *PetscOptionsObject)
412b3506946SBarry Smith {
413b3506946SBarry Smith   PetscErrorCode  ierr;
4144416b707SBarry Smith   PetscOptionItem next     = PetscOptionsObject->next;
415d5649816SBarry Smith   static int      mancount = 0;
416b3506946SBarry Smith   char            options[16];
4177aab2a10SBarry Smith   PetscBool       changedmethod = PETSC_FALSE;
418a23eb982SSurtai Han   PetscBool       stopasking    = PETSC_FALSE;
41988a9752cSBarry Smith   char            manname[16],textname[16];
4202657e9d9SBarry Smith   char            dir[1024];
421b3506946SBarry Smith 
4222a409bb0SBarry Smith   PetscFunctionBegin;
423b3506946SBarry Smith   /* the next line is a bug, this will only work if all processors are here, the comm passed in is ignored!!! */
424b3506946SBarry Smith   sprintf(options,"Options_%d",count++);
425a297a907SKarl Rupp 
4267325c4a4SBarry Smith   PetscOptionsObject->pprefix = PetscOptionsObject->prefix; /* SAWs will change this, so cannot pass prefix directly */
4271bc75a8dSBarry Smith 
428d50831c4SBarry Smith   ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s","_title");CHKERRQ(ierr);
42983355fc5SBarry Smith   PetscStackCallSAWs(SAWs_Register,(dir,&PetscOptionsObject->title,1,SAWs_READ,SAWs_STRING));
4307781c08eSBarry Smith   ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s","prefix");CHKERRQ(ierr);
43183355fc5SBarry Smith   PetscStackCallSAWs(SAWs_Register,(dir,&PetscOptionsObject->pprefix,1,SAWs_READ,SAWs_STRING));
4322657e9d9SBarry Smith   PetscStackCallSAWs(SAWs_Register,("/PETSc/Options/ChangedMethod",&changedmethod,1,SAWs_WRITE,SAWs_BOOLEAN));
433a23eb982SSurtai Han   PetscStackCallSAWs(SAWs_Register,("/PETSc/Options/StopAsking",&stopasking,1,SAWs_WRITE,SAWs_BOOLEAN));
434b3506946SBarry Smith 
435b3506946SBarry Smith   while (next) {
436d50831c4SBarry Smith     sprintf(manname,"_man_%d",mancount);
4372657e9d9SBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",manname);CHKERRQ(ierr);
4387781c08eSBarry Smith     PetscStackCallSAWs(SAWs_Register,(dir,&next->man,1,SAWs_READ,SAWs_STRING));
439d50831c4SBarry Smith     sprintf(textname,"_text_%d",mancount++);
4407781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",textname);CHKERRQ(ierr);
4417781c08eSBarry Smith     PetscStackCallSAWs(SAWs_Register,(dir,&next->text,1,SAWs_READ,SAWs_STRING));
4429f32e415SBarry Smith 
443b3506946SBarry Smith     switch (next->type) {
444b3506946SBarry Smith     case OPTION_HEAD:
445b3506946SBarry Smith       break;
446b3506946SBarry Smith     case OPTION_INT_ARRAY:
4477781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4482657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_INT));
449b3506946SBarry Smith       break;
450b3506946SBarry Smith     case OPTION_REAL_ARRAY:
4517781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4522657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_DOUBLE));
453b3506946SBarry Smith       break;
454b3506946SBarry Smith     case OPTION_INT:
4557781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4562657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_INT));
457b3506946SBarry Smith       break;
458b3506946SBarry Smith     case OPTION_REAL:
4597781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4602657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_DOUBLE));
461b3506946SBarry Smith       break;
4627781c08eSBarry Smith     case OPTION_BOOL:
4637781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4642657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_BOOLEAN));
4651ae3d29cSBarry Smith       break;
4667781c08eSBarry Smith     case OPTION_BOOL_ARRAY:
4677781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4682657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_BOOLEAN));
46971f08665SBarry Smith       break;
470b3506946SBarry Smith     case OPTION_STRING:
4717781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4727781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
4731ae3d29cSBarry Smith       break;
4741ae3d29cSBarry Smith     case OPTION_STRING_ARRAY:
4757781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4762657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_STRING));
477b3506946SBarry Smith       break;
478a264d7a6SBarry Smith     case OPTION_FLIST:
479a264d7a6SBarry Smith       {
480a264d7a6SBarry Smith       PetscInt ntext;
4817781c08eSBarry Smith       ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4827781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
483a264d7a6SBarry Smith       ierr = PetscFunctionListGet(next->flist,(const char***)&next->edata,&ntext);CHKERRQ(ierr);
484a264d7a6SBarry Smith       PetscStackCallSAWs(SAWs_Set_Legal_Variable_Values,(dir,ntext,next->edata));
485a264d7a6SBarry Smith       }
486a264d7a6SBarry Smith       break;
4871ae3d29cSBarry Smith     case OPTION_ELIST:
488a264d7a6SBarry Smith       {
489a264d7a6SBarry Smith       PetscInt ntext = next->nlist;
4907781c08eSBarry Smith       ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4917781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
492ead66b60SBarry Smith       ierr = PetscMalloc1((ntext+1),(char***)&next->edata);CHKERRQ(ierr);
4931ae3d29cSBarry Smith       ierr = PetscMemcpy(next->edata,next->list,ntext*sizeof(char*));CHKERRQ(ierr);
494a264d7a6SBarry Smith       PetscStackCallSAWs(SAWs_Set_Legal_Variable_Values,(dir,ntext,next->edata));
495a264d7a6SBarry Smith       }
496a264d7a6SBarry Smith       break;
497b3506946SBarry Smith     default:
498b3506946SBarry Smith       break;
499b3506946SBarry Smith     }
500b3506946SBarry Smith     next = next->next;
501b3506946SBarry Smith   }
502b3506946SBarry Smith 
503b3506946SBarry Smith   /* wait until accessor has unlocked the memory */
50464bbc9efSBarry Smith   PetscStackCallSAWs(SAWs_Push_Header,("index.html",OptionsHeader));
50564bbc9efSBarry Smith   PetscStackCallSAWs(SAWs_Push_Body,("index.html",2,OptionsBodyBottom));
5067aab2a10SBarry Smith   ierr = PetscSAWsBlock();CHKERRQ(ierr);
50764bbc9efSBarry Smith   PetscStackCallSAWs(SAWs_Pop_Header,("index.html"));
50864bbc9efSBarry Smith   PetscStackCallSAWs(SAWs_Pop_Body,("index.html",2));
509b3506946SBarry Smith 
51088a9752cSBarry Smith   /* determine if any values have been set in GUI */
51183355fc5SBarry Smith   next = PetscOptionsObject->next;
51288a9752cSBarry Smith   while (next) {
51388a9752cSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
514f7b25cf6SBarry Smith     PetscStackCallSAWs(SAWs_Selected,(dir,(int*)&next->set));
51588a9752cSBarry Smith     next = next->next;
51688a9752cSBarry Smith   }
51788a9752cSBarry Smith 
518b3506946SBarry Smith   /* reset counter to -2; this updates the screen with the new options for the selected method */
519f7b25cf6SBarry Smith   if (changedmethod) PetscOptionsObject->count = -2;
520b3506946SBarry Smith 
521a23eb982SSurtai Han   if (stopasking) {
522a23eb982SSurtai Han     PetscOptionsPublish      = PETSC_FALSE;
52312655325SBarry Smith     PetscOptionsObject->count = 0;//do not ask for same thing again
524a23eb982SSurtai Han   }
525a23eb982SSurtai Han 
5269a492a5cSBarry Smith   PetscStackCallSAWs(SAWs_Delete,("/PETSc/Options"));
527b3506946SBarry Smith   PetscFunctionReturn(0);
528b3506946SBarry Smith }
529b3506946SBarry Smith #endif
530b3506946SBarry Smith 
5316356e834SBarry Smith #undef __FUNCT__
53253acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEnd_Private"
5334416b707SBarry Smith PetscErrorCode PetscOptionsEnd_Private(PetscOptionItems *PetscOptionsObject)
53453acd3b1SBarry Smith {
53553acd3b1SBarry Smith   PetscErrorCode  ierr;
5364416b707SBarry Smith   PetscOptionItem last;
5376356e834SBarry Smith   char            option[256],value[1024],tmp[32];
538330cf3c9SBarry Smith   size_t          j;
53953acd3b1SBarry Smith 
54053acd3b1SBarry Smith   PetscFunctionBegin;
54183355fc5SBarry Smith   if (PetscOptionsObject->next) {
54283355fc5SBarry Smith     if (!PetscOptionsObject->count) {
543a264d7a6SBarry Smith #if defined(PETSC_HAVE_SAWS)
544f7b25cf6SBarry Smith       ierr = PetscOptionsSAWsInput(PetscOptionsObject);CHKERRQ(ierr);
545b3506946SBarry Smith #else
546e55864a3SBarry Smith       ierr = PetscOptionsGetFromTextInput(PetscOptionsObject);CHKERRQ(ierr);
547b3506946SBarry Smith #endif
548aee2cecaSBarry Smith     }
549aee2cecaSBarry Smith   }
5506356e834SBarry Smith 
551e55864a3SBarry Smith   ierr = PetscFree(PetscOptionsObject->title);CHKERRQ(ierr);
5526356e834SBarry Smith 
553e26ddf31SBarry Smith   /* reset counter to -2; this updates the screen with the new options for the selected method */
554e55864a3SBarry Smith   if (PetscOptionsObject->changedmethod) PetscOptionsObject->count = -2;
5557a72a596SBarry Smith   /* reset alreadyprinted flag */
556e55864a3SBarry Smith   PetscOptionsObject->alreadyprinted = PETSC_FALSE;
557e55864a3SBarry Smith   if (PetscOptionsObject->object) PetscOptionsObject->object->optionsprinted = PETSC_TRUE;
558e55864a3SBarry Smith   PetscOptionsObject->object = NULL;
55953acd3b1SBarry Smith 
560e55864a3SBarry Smith   while (PetscOptionsObject->next) {
561e55864a3SBarry Smith     if (PetscOptionsObject->next->set) {
562e55864a3SBarry Smith       if (PetscOptionsObject->prefix) {
56353acd3b1SBarry Smith         ierr = PetscStrcpy(option,"-");CHKERRQ(ierr);
564e55864a3SBarry Smith         ierr = PetscStrcat(option,PetscOptionsObject->prefix);CHKERRQ(ierr);
565e55864a3SBarry Smith         ierr = PetscStrcat(option,PetscOptionsObject->next->option+1);CHKERRQ(ierr);
5666356e834SBarry Smith       } else {
567e55864a3SBarry Smith         ierr = PetscStrcpy(option,PetscOptionsObject->next->option);CHKERRQ(ierr);
5686356e834SBarry Smith       }
5696356e834SBarry Smith 
570e55864a3SBarry Smith       switch (PetscOptionsObject->next->type) {
5716356e834SBarry Smith       case OPTION_HEAD:
5726356e834SBarry Smith         break;
5736356e834SBarry Smith       case OPTION_INT_ARRAY:
574e55864a3SBarry Smith         sprintf(value,"%d",(int)((PetscInt*)PetscOptionsObject->next->data)[0]);
575e55864a3SBarry Smith         for (j=1; j<PetscOptionsObject->next->arraylength; j++) {
576e55864a3SBarry Smith           sprintf(tmp,"%d",(int)((PetscInt*)PetscOptionsObject->next->data)[j]);
5776356e834SBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5786356e834SBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5796356e834SBarry Smith         }
5806356e834SBarry Smith         break;
5816356e834SBarry Smith       case OPTION_INT:
582e55864a3SBarry Smith         sprintf(value,"%d",(int) *(PetscInt*)PetscOptionsObject->next->data);
5836356e834SBarry Smith         break;
5846356e834SBarry Smith       case OPTION_REAL:
585e55864a3SBarry Smith         sprintf(value,"%g",(double) *(PetscReal*)PetscOptionsObject->next->data);
5866356e834SBarry Smith         break;
5876356e834SBarry Smith       case OPTION_REAL_ARRAY:
588e55864a3SBarry Smith         sprintf(value,"%g",(double)((PetscReal*)PetscOptionsObject->next->data)[0]);
589e55864a3SBarry Smith         for (j=1; j<PetscOptionsObject->next->arraylength; j++) {
590e55864a3SBarry Smith           sprintf(tmp,"%g",(double)((PetscReal*)PetscOptionsObject->next->data)[j]);
5916356e834SBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5926356e834SBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5936356e834SBarry Smith         }
5946356e834SBarry Smith         break;
595050cccc3SHong Zhang       case OPTION_SCALAR_ARRAY:
59695f3a755SJose E. Roman         sprintf(value,"%g+%gi",(double)PetscRealPart(((PetscScalar*)PetscOptionsObject->next->data)[0]),(double)PetscImaginaryPart(((PetscScalar*)PetscOptionsObject->next->data)[0]));
597050cccc3SHong Zhang         for (j=1; j<PetscOptionsObject->next->arraylength; j++) {
59895f3a755SJose E. Roman           sprintf(tmp,"%g+%gi",(double)PetscRealPart(((PetscScalar*)PetscOptionsObject->next->data)[j]),(double)PetscImaginaryPart(((PetscScalar*)PetscOptionsObject->next->data)[j]));
599050cccc3SHong Zhang           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
600050cccc3SHong Zhang           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
601050cccc3SHong Zhang         }
602050cccc3SHong Zhang         break;
6037781c08eSBarry Smith       case OPTION_BOOL:
604e55864a3SBarry Smith         sprintf(value,"%d",*(int*)PetscOptionsObject->next->data);
6056356e834SBarry Smith         break;
6067781c08eSBarry Smith       case OPTION_BOOL_ARRAY:
607e55864a3SBarry Smith         sprintf(value,"%d",(int)((PetscBool*)PetscOptionsObject->next->data)[0]);
608e55864a3SBarry Smith         for (j=1; j<PetscOptionsObject->next->arraylength; j++) {
609e55864a3SBarry Smith           sprintf(tmp,"%d",(int)((PetscBool*)PetscOptionsObject->next->data)[j]);
6101ae3d29cSBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
6111ae3d29cSBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
6121ae3d29cSBarry Smith         }
6131ae3d29cSBarry Smith         break;
614a264d7a6SBarry Smith       case OPTION_FLIST:
6156991f827SBarry Smith         ierr = PetscStrcpy(value,(char*)PetscOptionsObject->next->data);CHKERRQ(ierr);
6166991f827SBarry Smith         break;
6171ae3d29cSBarry Smith       case OPTION_ELIST:
618e55864a3SBarry Smith         ierr = PetscStrcpy(value,(char*)PetscOptionsObject->next->data);CHKERRQ(ierr);
6196356e834SBarry Smith         break;
6201ae3d29cSBarry Smith       case OPTION_STRING:
621e55864a3SBarry Smith         ierr = PetscStrcpy(value,(char*)PetscOptionsObject->next->data);CHKERRQ(ierr);
622d51da6bfSBarry Smith         break;
6231ae3d29cSBarry Smith       case OPTION_STRING_ARRAY:
624e55864a3SBarry Smith         sprintf(value,"%s",((char**)PetscOptionsObject->next->data)[0]);
625e55864a3SBarry Smith         for (j=1; j<PetscOptionsObject->next->arraylength; j++) {
626e55864a3SBarry Smith           sprintf(tmp,"%s",((char**)PetscOptionsObject->next->data)[j]);
6271ae3d29cSBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
6281ae3d29cSBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
6291ae3d29cSBarry Smith         }
6306356e834SBarry Smith         break;
6316356e834SBarry Smith       }
632c5929fdfSBarry Smith       ierr = PetscOptionsSetValue(PetscOptionsObject->options,option,value);CHKERRQ(ierr);
6336356e834SBarry Smith     }
6346991f827SBarry Smith     if (PetscOptionsObject->next->type == OPTION_ELIST) {
6356991f827SBarry Smith       ierr = PetscStrNArrayDestroy(PetscOptionsObject->next->nlist,(char ***)&PetscOptionsObject->next->list);CHKERRQ(ierr);
6366991f827SBarry Smith     }
637e55864a3SBarry Smith     ierr   = PetscFree(PetscOptionsObject->next->text);CHKERRQ(ierr);
638e55864a3SBarry Smith     ierr   = PetscFree(PetscOptionsObject->next->option);CHKERRQ(ierr);
639e55864a3SBarry Smith     ierr   = PetscFree(PetscOptionsObject->next->man);CHKERRQ(ierr);
640e55864a3SBarry Smith     ierr   = PetscFree(PetscOptionsObject->next->edata);CHKERRQ(ierr);
641c979a496SBarry Smith 
64283355fc5SBarry Smith     if ((PetscOptionsObject->next->type == OPTION_STRING) || (PetscOptionsObject->next->type == OPTION_FLIST) || (PetscOptionsObject->next->type == OPTION_ELIST)){
64383355fc5SBarry Smith       free(PetscOptionsObject->next->data);
644c979a496SBarry Smith     } else {
64583355fc5SBarry Smith       ierr   = PetscFree(PetscOptionsObject->next->data);CHKERRQ(ierr);
646c979a496SBarry Smith     }
6477781c08eSBarry Smith 
64883355fc5SBarry Smith     last                    = PetscOptionsObject->next;
64983355fc5SBarry Smith     PetscOptionsObject->next = PetscOptionsObject->next->next;
6506356e834SBarry Smith     ierr                    = PetscFree(last);CHKERRQ(ierr);
6516356e834SBarry Smith   }
652f59f755dSBarry Smith   ierr = PetscFree(PetscOptionsObject->prefix);CHKERRQ(ierr);
653e55864a3SBarry Smith   PetscOptionsObject->next = 0;
65453acd3b1SBarry Smith   PetscFunctionReturn(0);
65553acd3b1SBarry Smith }
65653acd3b1SBarry Smith 
65753acd3b1SBarry Smith #undef __FUNCT__
658e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsEnum_Private"
65953acd3b1SBarry Smith /*@C
66053acd3b1SBarry Smith    PetscOptionsEnum - Gets the enum value for a particular option in the database.
66153acd3b1SBarry Smith 
6623f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
66353acd3b1SBarry Smith 
66453acd3b1SBarry Smith    Input Parameters:
66553acd3b1SBarry Smith +  opt - option name
66653acd3b1SBarry Smith .  text - short string that describes the option
66753acd3b1SBarry Smith .  man - manual page with additional information on option
66853acd3b1SBarry Smith .  list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null
6690fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
6700fdccdaeSBarry Smith $                 PetscOptionsEnum(..., obj->value,&object->value,...) or
6710fdccdaeSBarry Smith $                 value = defaultvalue
6720fdccdaeSBarry Smith $                 PetscOptionsEnum(..., value,&value,&flg);
6730fdccdaeSBarry Smith $                 if (flg) {
67453acd3b1SBarry Smith 
67553acd3b1SBarry Smith    Output Parameter:
67653acd3b1SBarry Smith +  value - the  value to return
677b32e0204SMatthew G Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
67853acd3b1SBarry Smith 
67953acd3b1SBarry Smith    Level: beginner
68053acd3b1SBarry Smith 
68153acd3b1SBarry Smith    Concepts: options database
68253acd3b1SBarry Smith 
68353acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
68453acd3b1SBarry Smith 
68553acd3b1SBarry Smith           list is usually something like PCASMTypes or some other predefined list of enum names
68653acd3b1SBarry Smith 
68753acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
688acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
689acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
69053acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
69153acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
692acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
693a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
69453acd3b1SBarry Smith @*/
6954416b707SBarry 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)
69653acd3b1SBarry Smith {
69753acd3b1SBarry Smith   PetscErrorCode ierr;
69853acd3b1SBarry Smith   PetscInt       ntext = 0;
699aa5bb8c0SSatish Balay   PetscInt       tval;
700ace3abfcSBarry Smith   PetscBool      tflg;
70153acd3b1SBarry Smith 
70253acd3b1SBarry Smith   PetscFunctionBegin;
70353acd3b1SBarry Smith   while (list[ntext++]) {
704e32f2f54SBarry Smith     if (ntext > 50) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument appears to be wrong or have more than 50 entries");
70553acd3b1SBarry Smith   }
706e32f2f54SBarry Smith   if (ntext < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument must have at least two entries: typename and type prefix");
70753acd3b1SBarry Smith   ntext -= 3;
708e55864a3SBarry Smith   ierr   = PetscOptionsEList_Private(PetscOptionsObject,opt,text,man,list,ntext,list[currentvalue],&tval,&tflg);CHKERRQ(ierr);
709aa5bb8c0SSatish Balay   /* with PETSC_USE_64BIT_INDICES sizeof(PetscInt) != sizeof(PetscEnum) */
710aa5bb8c0SSatish Balay   if (tflg) *value = (PetscEnum)tval;
711aa5bb8c0SSatish Balay   if (set)  *set   = tflg;
71253acd3b1SBarry Smith   PetscFunctionReturn(0);
71353acd3b1SBarry Smith }
71453acd3b1SBarry Smith 
715d3e47460SLisandro Dalcin #undef __FUNCT__
716d3e47460SLisandro Dalcin #define __FUNCT__ "PetscOptionsEnumArray_Private"
717d3e47460SLisandro Dalcin /*@C
718d3e47460SLisandro Dalcin    PetscOptionsEnumArray - Gets an array of enum values for a particular
719d3e47460SLisandro Dalcin    option in the database.
720d3e47460SLisandro Dalcin 
721d3e47460SLisandro Dalcin    Logically Collective on the communicator passed in PetscOptionsBegin()
722d3e47460SLisandro Dalcin 
723d3e47460SLisandro Dalcin    Input Parameters:
724d3e47460SLisandro Dalcin +  opt - the option one is seeking
725d3e47460SLisandro Dalcin .  text - short string describing option
726d3e47460SLisandro Dalcin .  man - manual page for option
727d3e47460SLisandro Dalcin -  n - maximum number of values
728d3e47460SLisandro Dalcin 
729d3e47460SLisandro Dalcin    Output Parameter:
730d3e47460SLisandro Dalcin +  value - location to copy values
731d3e47460SLisandro Dalcin .  n - actual number of values found
732d3e47460SLisandro Dalcin -  set - PETSC_TRUE if found, else PETSC_FALSE
733d3e47460SLisandro Dalcin 
734d3e47460SLisandro Dalcin    Level: beginner
735d3e47460SLisandro Dalcin 
736d3e47460SLisandro Dalcin    Notes:
737d3e47460SLisandro Dalcin    The array must be passed as a comma separated list.
738d3e47460SLisandro Dalcin 
739d3e47460SLisandro Dalcin    There must be no intervening spaces between the values.
740d3e47460SLisandro Dalcin 
741d3e47460SLisandro Dalcin    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
742d3e47460SLisandro Dalcin 
743d3e47460SLisandro Dalcin    Concepts: options database^array of enums
744d3e47460SLisandro Dalcin 
745d3e47460SLisandro Dalcin .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
746d3e47460SLisandro Dalcin           PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
747d3e47460SLisandro Dalcin           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
748d3e47460SLisandro Dalcin           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
749d3e47460SLisandro Dalcin           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
750d3e47460SLisandro Dalcin           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsRealArray()
751d3e47460SLisandro Dalcin @*/
7524416b707SBarry 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)
753d3e47460SLisandro Dalcin {
754d3e47460SLisandro Dalcin   PetscInt        i,nlist = 0;
7554416b707SBarry Smith   PetscOptionItem amsopt;
756d3e47460SLisandro Dalcin   PetscErrorCode  ierr;
757d3e47460SLisandro Dalcin 
758d3e47460SLisandro Dalcin   PetscFunctionBegin;
759d3e47460SLisandro 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");
760d3e47460SLisandro Dalcin   if (nlist < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument must have at least two entries: typename and type prefix");
761d3e47460SLisandro Dalcin   nlist -= 3; /* drop enum name, prefix, and null termination */
762d3e47460SLisandro Dalcin   if (0 && !PetscOptionsObject->count) { /* XXX Requires additional support */
763d3e47460SLisandro Dalcin     PetscEnum *vals;
7644416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_INT_ARRAY/*XXX OPTION_ENUM_ARRAY*/,&amsopt);CHKERRQ(ierr);
765d3e47460SLisandro Dalcin     ierr = PetscStrNArrayallocpy(nlist,list,(char***)&amsopt->list);CHKERRQ(ierr);
766d3e47460SLisandro Dalcin     amsopt->nlist = nlist;
767d3e47460SLisandro Dalcin     ierr = PetscMalloc1(*n,(PetscEnum**)&amsopt->data);CHKERRQ(ierr);
768d3e47460SLisandro Dalcin     amsopt->arraylength = *n;
769d3e47460SLisandro Dalcin     vals = (PetscEnum*)amsopt->data;
770d3e47460SLisandro Dalcin     for (i=0; i<*n; i++) vals[i] = value[i];
771d3e47460SLisandro Dalcin   }
772c5929fdfSBarry Smith   ierr = PetscOptionsGetEnumArray(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,list,value,n,set);CHKERRQ(ierr);
773d3e47460SLisandro Dalcin   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
774d3e47460SLisandro Dalcin     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%s",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,list[value[0]]);CHKERRQ(ierr);
775d3e47460SLisandro Dalcin     for (i=1; i<*n; i++) {ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%s",list[value[i]]);CHKERRQ(ierr);}
776d3e47460SLisandro Dalcin     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (choose from)",text);CHKERRQ(ierr);
777d3e47460SLisandro Dalcin     for (i=0; i<nlist; i++) {ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," %s",list[i]);CHKERRQ(ierr);}
778d3e47460SLisandro Dalcin     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," (%s)\n",ManSection(man));CHKERRQ(ierr);
779d3e47460SLisandro Dalcin   }
780d3e47460SLisandro Dalcin   PetscFunctionReturn(0);
781d3e47460SLisandro Dalcin }
782d3e47460SLisandro Dalcin 
78353acd3b1SBarry Smith /* -------------------------------------------------------------------------------------------------------------*/
78453acd3b1SBarry Smith #undef __FUNCT__
785e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsInt_Private"
78653acd3b1SBarry Smith /*@C
78753acd3b1SBarry Smith    PetscOptionsInt - Gets the integer value for a particular option in the database.
78853acd3b1SBarry Smith 
7893f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
79053acd3b1SBarry Smith 
79153acd3b1SBarry Smith    Input Parameters:
79253acd3b1SBarry Smith +  opt - option name
79353acd3b1SBarry Smith .  text - short string that describes the option
79453acd3b1SBarry Smith .  man - manual page with additional information on option
7950fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
7960fdccdaeSBarry Smith $                 PetscOptionsInt(..., obj->value,&object->value,...) or
7970fdccdaeSBarry Smith $                 value = defaultvalue
7980fdccdaeSBarry Smith $                 PetscOptionsInt(..., value,&value,&flg);
7990fdccdaeSBarry Smith $                 if (flg) {
80053acd3b1SBarry Smith 
80153acd3b1SBarry Smith    Output Parameter:
80253acd3b1SBarry Smith +  value - the integer value to return
80353acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
80453acd3b1SBarry Smith 
80553acd3b1SBarry Smith    Level: beginner
80653acd3b1SBarry Smith 
80753acd3b1SBarry Smith    Concepts: options database^has int
80853acd3b1SBarry Smith 
80953acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
81053acd3b1SBarry Smith 
81153acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
812acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
813acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
81453acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
81553acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
816acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
817a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
81853acd3b1SBarry Smith @*/
8194416b707SBarry Smith PetscErrorCode  PetscOptionsInt_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscInt currentvalue,PetscInt *value,PetscBool  *set)
82053acd3b1SBarry Smith {
82153acd3b1SBarry Smith   PetscErrorCode  ierr;
8224416b707SBarry Smith   PetscOptionItem amsopt;
82312655325SBarry Smith   PetscBool       wasset;
82453acd3b1SBarry Smith 
82553acd3b1SBarry Smith   PetscFunctionBegin;
826e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
8274416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_INT,&amsopt);CHKERRQ(ierr);
8286356e834SBarry Smith     ierr = PetscMalloc(sizeof(PetscInt),&amsopt->data);CHKERRQ(ierr);
82912655325SBarry Smith     *(PetscInt*)amsopt->data = currentvalue;
8303e211508SBarry Smith 
831c5929fdfSBarry Smith     ierr = PetscOptionsGetInt(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,&currentvalue,&wasset);CHKERRQ(ierr);
8323e211508SBarry Smith     if (wasset) {
83312655325SBarry Smith       *(PetscInt*)amsopt->data = currentvalue;
8343e211508SBarry Smith     }
835af6d86caSBarry Smith   }
836c5929fdfSBarry Smith   ierr = PetscOptionsGetInt(NULL,PetscOptionsObject->prefix,opt,value,set);CHKERRQ(ierr);
837e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
8381a1499c8SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%d>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,currentvalue,text,ManSection(man));CHKERRQ(ierr);
83953acd3b1SBarry Smith   }
84053acd3b1SBarry Smith   PetscFunctionReturn(0);
84153acd3b1SBarry Smith }
84253acd3b1SBarry Smith 
84353acd3b1SBarry Smith #undef __FUNCT__
844e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsString_Private"
84553acd3b1SBarry Smith /*@C
84653acd3b1SBarry Smith    PetscOptionsString - Gets the string value for a particular option in the database.
84753acd3b1SBarry Smith 
8483f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
84953acd3b1SBarry Smith 
85053acd3b1SBarry Smith    Input Parameters:
85153acd3b1SBarry Smith +  opt - option name
85253acd3b1SBarry Smith .  text - short string that describes the option
85353acd3b1SBarry Smith .  man - manual page with additional information on option
8540fdccdaeSBarry Smith .  currentvalue - the current value; caller is responsible for setting this value correctly. This is not used to set value
855bcbf2dc5SJed Brown -  len - length of the result string including null terminator
85653acd3b1SBarry Smith 
85753acd3b1SBarry Smith    Output Parameter:
85853acd3b1SBarry Smith +  value - the value to return
85953acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
86053acd3b1SBarry Smith 
86153acd3b1SBarry Smith    Level: beginner
86253acd3b1SBarry Smith 
86353acd3b1SBarry Smith    Concepts: options database^has int
86453acd3b1SBarry Smith 
86553acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
86653acd3b1SBarry Smith 
8677fccdfe4SBarry 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).
8687fccdfe4SBarry Smith 
869c5929fdfSBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(NULL,),
870acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
871acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsReal(), PetscOptionsBool(),
87253acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
87353acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
874acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
875a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
87653acd3b1SBarry Smith @*/
8774416b707SBarry 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)
87853acd3b1SBarry Smith {
87953acd3b1SBarry Smith   PetscErrorCode  ierr;
8804416b707SBarry Smith   PetscOptionItem amsopt;
88153acd3b1SBarry Smith 
88253acd3b1SBarry Smith   PetscFunctionBegin;
8831a1499c8SBarry Smith   if (!PetscOptionsObject->count) {
8844416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
88564facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
8860fdccdaeSBarry Smith     ierr = PetscStrdup(currentvalue ? currentvalue : "",(char**)&amsopt->data);CHKERRQ(ierr);
887af6d86caSBarry Smith   }
888c5929fdfSBarry Smith   ierr = PetscOptionsGetString(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,len,set);CHKERRQ(ierr);
889e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
8901a1499c8SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%s>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,currentvalue,text,ManSection(man));CHKERRQ(ierr);
89153acd3b1SBarry Smith   }
89253acd3b1SBarry Smith   PetscFunctionReturn(0);
89353acd3b1SBarry Smith }
89453acd3b1SBarry Smith 
89553acd3b1SBarry Smith #undef __FUNCT__
896e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsReal_Private"
89753acd3b1SBarry Smith /*@C
89853acd3b1SBarry Smith    PetscOptionsReal - Gets the PetscReal value for a particular option in the database.
89953acd3b1SBarry Smith 
9003f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
90153acd3b1SBarry Smith 
90253acd3b1SBarry Smith    Input Parameters:
90353acd3b1SBarry Smith +  opt - option name
90453acd3b1SBarry Smith .  text - short string that describes the option
90553acd3b1SBarry Smith .  man - manual page with additional information on option
9060fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
9070fdccdaeSBarry Smith $                 PetscOptionsReal(..., obj->value,&object->value,...) or
9080fdccdaeSBarry Smith $                 value = defaultvalue
9090fdccdaeSBarry Smith $                 PetscOptionsReal(..., value,&value,&flg);
9100fdccdaeSBarry Smith $                 if (flg) {
91153acd3b1SBarry Smith 
91253acd3b1SBarry Smith    Output Parameter:
91353acd3b1SBarry Smith +  value - the value to return
91453acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
91553acd3b1SBarry Smith 
91653acd3b1SBarry Smith    Level: beginner
91753acd3b1SBarry Smith 
91853acd3b1SBarry Smith    Concepts: options database^has int
91953acd3b1SBarry Smith 
92053acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
92153acd3b1SBarry Smith 
922c5929fdfSBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(NULL,),
923acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
924acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
92553acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
92653acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
927acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
928a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
92953acd3b1SBarry Smith @*/
9304416b707SBarry Smith PetscErrorCode  PetscOptionsReal_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscReal currentvalue,PetscReal *value,PetscBool  *set)
93153acd3b1SBarry Smith {
93253acd3b1SBarry Smith   PetscErrorCode  ierr;
9334416b707SBarry Smith   PetscOptionItem amsopt;
93453acd3b1SBarry Smith 
93553acd3b1SBarry Smith   PetscFunctionBegin;
936e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
9374416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_REAL,&amsopt);CHKERRQ(ierr);
938538aa990SBarry Smith     ierr = PetscMalloc(sizeof(PetscReal),&amsopt->data);CHKERRQ(ierr);
939a297a907SKarl Rupp 
9400fdccdaeSBarry Smith     *(PetscReal*)amsopt->data = currentvalue;
941538aa990SBarry Smith   }
942c5929fdfSBarry Smith   ierr = PetscOptionsGetReal(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,set);CHKERRQ(ierr);
9431a1499c8SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
9441a1499c8SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%g>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,(double)currentvalue,text,ManSection(man));CHKERRQ(ierr);
94553acd3b1SBarry Smith   }
94653acd3b1SBarry Smith   PetscFunctionReturn(0);
94753acd3b1SBarry Smith }
94853acd3b1SBarry Smith 
94953acd3b1SBarry Smith #undef __FUNCT__
950e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsScalar_Private"
95153acd3b1SBarry Smith /*@C
95253acd3b1SBarry Smith    PetscOptionsScalar - Gets the scalar value for a particular option in the database.
95353acd3b1SBarry Smith 
9543f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
95553acd3b1SBarry Smith 
95653acd3b1SBarry Smith    Input Parameters:
95753acd3b1SBarry Smith +  opt - option name
95853acd3b1SBarry Smith .  text - short string that describes the option
95953acd3b1SBarry Smith .  man - manual page with additional information on option
9600fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
9610fdccdaeSBarry Smith $                 PetscOptionsScalar(..., obj->value,&object->value,...) or
9620fdccdaeSBarry Smith $                 value = defaultvalue
9630fdccdaeSBarry Smith $                 PetscOptionsScalar(..., value,&value,&flg);
9640fdccdaeSBarry Smith $                 if (flg) {
9650fdccdaeSBarry Smith 
96653acd3b1SBarry Smith 
96753acd3b1SBarry Smith    Output Parameter:
96853acd3b1SBarry Smith +  value - the value to return
96953acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
97053acd3b1SBarry Smith 
97153acd3b1SBarry Smith    Level: beginner
97253acd3b1SBarry Smith 
97353acd3b1SBarry Smith    Concepts: options database^has int
97453acd3b1SBarry Smith 
97553acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
97653acd3b1SBarry Smith 
977c5929fdfSBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(NULL,),
978acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
979acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
98053acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
98153acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
982acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
983a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
98453acd3b1SBarry Smith @*/
9854416b707SBarry Smith PetscErrorCode  PetscOptionsScalar_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscScalar currentvalue,PetscScalar *value,PetscBool  *set)
98653acd3b1SBarry Smith {
98753acd3b1SBarry Smith   PetscErrorCode ierr;
98853acd3b1SBarry Smith 
98953acd3b1SBarry Smith   PetscFunctionBegin;
99053acd3b1SBarry Smith #if !defined(PETSC_USE_COMPLEX)
9910fdccdaeSBarry Smith   ierr = PetscOptionsReal(opt,text,man,currentvalue,value,set);CHKERRQ(ierr);
99253acd3b1SBarry Smith #else
993c5929fdfSBarry Smith   ierr = PetscOptionsGetScalar(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,set);CHKERRQ(ierr);
99453acd3b1SBarry Smith #endif
99553acd3b1SBarry Smith   PetscFunctionReturn(0);
99653acd3b1SBarry Smith }
99753acd3b1SBarry Smith 
99853acd3b1SBarry Smith #undef __FUNCT__
999e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsName_Private"
100053acd3b1SBarry Smith /*@C
100190d69ab7SBarry 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
100290d69ab7SBarry Smith                       its value is set to false.
100353acd3b1SBarry Smith 
10043f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
100553acd3b1SBarry Smith 
100653acd3b1SBarry Smith    Input Parameters:
100753acd3b1SBarry Smith +  opt - option name
100853acd3b1SBarry Smith .  text - short string that describes the option
100953acd3b1SBarry Smith -  man - manual page with additional information on option
101053acd3b1SBarry Smith 
101153acd3b1SBarry Smith    Output Parameter:
101253acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
101353acd3b1SBarry Smith 
101453acd3b1SBarry Smith    Level: beginner
101553acd3b1SBarry Smith 
101653acd3b1SBarry Smith    Concepts: options database^has int
101753acd3b1SBarry Smith 
101853acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
101953acd3b1SBarry Smith 
1020c5929fdfSBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(NULL,),
1021acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
1022acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
102353acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
102453acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1025acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1026a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
102753acd3b1SBarry Smith @*/
10284416b707SBarry Smith PetscErrorCode  PetscOptionsName_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool  *flg)
102953acd3b1SBarry Smith {
103053acd3b1SBarry Smith   PetscErrorCode  ierr;
10314416b707SBarry Smith   PetscOptionItem amsopt;
103253acd3b1SBarry Smith 
103353acd3b1SBarry Smith   PetscFunctionBegin;
1034e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
10354416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1036ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1037a297a907SKarl Rupp 
1038ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
10391ae3d29cSBarry Smith   }
1040c5929fdfSBarry Smith   ierr = PetscOptionsHasName(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,flg);CHKERRQ(ierr);
1041e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1042e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
104353acd3b1SBarry Smith   }
104453acd3b1SBarry Smith   PetscFunctionReturn(0);
104553acd3b1SBarry Smith }
104653acd3b1SBarry Smith 
104753acd3b1SBarry Smith #undef __FUNCT__
1048e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsFList_Private"
104953acd3b1SBarry Smith /*@C
1050a264d7a6SBarry Smith      PetscOptionsFList - Puts a list of option values that a single one may be selected from
105153acd3b1SBarry Smith 
10523f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
105353acd3b1SBarry Smith 
105453acd3b1SBarry Smith    Input Parameters:
105553acd3b1SBarry Smith +  opt - option name
105653acd3b1SBarry Smith .  text - short string that describes the option
105753acd3b1SBarry Smith .  man - manual page with additional information on option
105853acd3b1SBarry Smith .  list - the possible choices
10590fdccdaeSBarry Smith .  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with
10600fdccdaeSBarry Smith $                 PetscOptionsFlist(..., obj->value,value,len,&flg);
10610fdccdaeSBarry Smith $                 if (flg) {
10623cc1e11dSBarry Smith -  len - the length of the character array value
106353acd3b1SBarry Smith 
106453acd3b1SBarry Smith    Output Parameter:
106553acd3b1SBarry Smith +  value - the value to return
106653acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
106753acd3b1SBarry Smith 
106853acd3b1SBarry Smith    Level: intermediate
106953acd3b1SBarry Smith 
107053acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
107153acd3b1SBarry Smith 
107253acd3b1SBarry Smith    See PetscOptionsEList() for when the choices are given in a string array
107353acd3b1SBarry Smith 
107453acd3b1SBarry Smith    To get a listing of all currently specified options,
107588c29154SBarry Smith     see PetscOptionsView() or PetscOptionsGetAll()
107653acd3b1SBarry Smith 
1077eabe10d7SBarry Smith    Developer Note: This cannot check for invalid selection because of things like MATAIJ that are not included in the list
1078eabe10d7SBarry Smith 
107953acd3b1SBarry Smith    Concepts: options database^list
108053acd3b1SBarry Smith 
1081c5929fdfSBarry Smith .seealso: PetscOptionsGetInt(NULL,), PetscOptionsGetReal(),
1082acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
108353acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
108453acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1085acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1086a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsEnum()
108753acd3b1SBarry Smith @*/
10884416b707SBarry 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)
108953acd3b1SBarry Smith {
109053acd3b1SBarry Smith   PetscErrorCode  ierr;
10914416b707SBarry Smith   PetscOptionItem amsopt;
109253acd3b1SBarry Smith 
109353acd3b1SBarry Smith   PetscFunctionBegin;
10941a1499c8SBarry Smith   if (!PetscOptionsObject->count) {
10954416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,ltext,man,OPTION_FLIST,&amsopt);CHKERRQ(ierr);
109664facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
10970fdccdaeSBarry Smith     ierr = PetscStrdup(currentvalue ? currentvalue : "",(char**)&amsopt->data);CHKERRQ(ierr);
10983cc1e11dSBarry Smith     amsopt->flist = list;
10993cc1e11dSBarry Smith   }
1100c5929fdfSBarry Smith   ierr = PetscOptionsGetString(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,len,set);CHKERRQ(ierr);
11011a1499c8SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
11021a1499c8SBarry Smith     ierr = PetscFunctionListPrintTypes(PetscOptionsObject->comm,stdout,PetscOptionsObject->prefix,opt,ltext,man,list,currentvalue);CHKERRQ(ierr);CHKERRQ(ierr);
110353acd3b1SBarry Smith   }
110453acd3b1SBarry Smith   PetscFunctionReturn(0);
110553acd3b1SBarry Smith }
110653acd3b1SBarry Smith 
110753acd3b1SBarry Smith #undef __FUNCT__
1108e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsEList_Private"
110953acd3b1SBarry Smith /*@C
111053acd3b1SBarry Smith      PetscOptionsEList - Puts a list of option values that a single one may be selected from
111153acd3b1SBarry Smith 
11123f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
111353acd3b1SBarry Smith 
111453acd3b1SBarry Smith    Input Parameters:
111553acd3b1SBarry Smith +  opt - option name
111653acd3b1SBarry Smith .  ltext - short string that describes the option
111753acd3b1SBarry Smith .  man - manual page with additional information on option
1118a264d7a6SBarry Smith .  list - the possible choices (one of these must be selected, anything else is invalid)
111953acd3b1SBarry Smith .  ntext - number of choices
11200fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with
11210fdccdaeSBarry Smith $                 PetscOptionsElist(..., obj->value,&value,&flg);
11220fdccdaeSBarry Smith $                 if (flg) {
11230fdccdaeSBarry Smith 
112453acd3b1SBarry Smith 
112553acd3b1SBarry Smith    Output Parameter:
112653acd3b1SBarry Smith +  value - the index of the value to return
112753acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
112853acd3b1SBarry Smith 
112953acd3b1SBarry Smith    Level: intermediate
113053acd3b1SBarry Smith 
113153acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
113253acd3b1SBarry Smith 
1133a264d7a6SBarry Smith    See PetscOptionsFList() for when the choices are given in a PetscFunctionList()
113453acd3b1SBarry Smith 
113553acd3b1SBarry Smith    Concepts: options database^list
113653acd3b1SBarry Smith 
1137c5929fdfSBarry Smith .seealso: PetscOptionsGetInt(NULL,), PetscOptionsGetReal(),
1138acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
113953acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
114053acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1141acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1142a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEnum()
114353acd3b1SBarry Smith @*/
11444416b707SBarry 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)
114553acd3b1SBarry Smith {
114653acd3b1SBarry Smith   PetscErrorCode  ierr;
114753acd3b1SBarry Smith   PetscInt        i;
11484416b707SBarry Smith   PetscOptionItem amsopt;
114953acd3b1SBarry Smith 
115053acd3b1SBarry Smith   PetscFunctionBegin;
11511a1499c8SBarry Smith   if (!PetscOptionsObject->count) {
11524416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,ltext,man,OPTION_ELIST,&amsopt);CHKERRQ(ierr);
115364facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
11540fdccdaeSBarry Smith     ierr = PetscStrdup(currentvalue ? currentvalue : "",(char**)&amsopt->data);CHKERRQ(ierr);
11556991f827SBarry Smith     ierr = PetscStrNArrayallocpy(ntext,list,(char***)&amsopt->list);CHKERRQ(ierr);
11561ae3d29cSBarry Smith     amsopt->nlist = ntext;
11571ae3d29cSBarry Smith   }
1158c5929fdfSBarry Smith   ierr = PetscOptionsGetEList(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,list,ntext,value,set);CHKERRQ(ierr);
11591a1499c8SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1160e3f729a5SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%s> %s (choose one of)",PetscOptionsObject->prefix?PetscOptionsObject->prefix:"",opt+1,currentvalue,ltext);CHKERRQ(ierr);
116153acd3b1SBarry Smith     for (i=0; i<ntext; i++) {
1162e55864a3SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," %s",list[i]);CHKERRQ(ierr);
116353acd3b1SBarry Smith     }
1164e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," (%s)\n",ManSection(man));CHKERRQ(ierr);
116553acd3b1SBarry Smith   }
116653acd3b1SBarry Smith   PetscFunctionReturn(0);
116753acd3b1SBarry Smith }
116853acd3b1SBarry Smith 
116953acd3b1SBarry Smith #undef __FUNCT__
1170e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsBoolGroupBegin_Private"
117153acd3b1SBarry Smith /*@C
1172acfcf0e5SJed Brown      PetscOptionsBoolGroupBegin - First in a series of logical queries on the options database for
1173d5649816SBarry Smith        which at most a single value can be true.
117453acd3b1SBarry Smith 
11753f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
117653acd3b1SBarry Smith 
117753acd3b1SBarry Smith    Input Parameters:
117853acd3b1SBarry Smith +  opt - option name
117953acd3b1SBarry Smith .  text - short string that describes the option
118053acd3b1SBarry Smith -  man - manual page with additional information on option
118153acd3b1SBarry Smith 
118253acd3b1SBarry Smith    Output Parameter:
118353acd3b1SBarry Smith .  flg - whether that option was set or not
118453acd3b1SBarry Smith 
118553acd3b1SBarry Smith    Level: intermediate
118653acd3b1SBarry Smith 
118753acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
118853acd3b1SBarry Smith 
1189acfcf0e5SJed Brown    Must be followed by 0 or more PetscOptionsBoolGroup()s and PetscOptionsBoolGroupEnd()
119053acd3b1SBarry Smith 
119153acd3b1SBarry Smith     Concepts: options database^logical group
119253acd3b1SBarry Smith 
1193c5929fdfSBarry Smith .seealso: PetscOptionsGetInt(NULL,), PetscOptionsGetReal(),
1194acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
119553acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
119653acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1197acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1198a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
119953acd3b1SBarry Smith @*/
12004416b707SBarry Smith PetscErrorCode  PetscOptionsBoolGroupBegin_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool  *flg)
120153acd3b1SBarry Smith {
120253acd3b1SBarry Smith   PetscErrorCode  ierr;
12034416b707SBarry Smith   PetscOptionItem amsopt;
120453acd3b1SBarry Smith 
120553acd3b1SBarry Smith   PetscFunctionBegin;
1206e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
12074416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1208ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1209a297a907SKarl Rupp 
1210ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
12111ae3d29cSBarry Smith   }
121268b16fdaSBarry Smith   *flg = PETSC_FALSE;
1213c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,flg,NULL);CHKERRQ(ierr);
1214e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1215e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  Pick at most one of -------------\n");CHKERRQ(ierr);
1216e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"    -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
121753acd3b1SBarry Smith   }
121853acd3b1SBarry Smith   PetscFunctionReturn(0);
121953acd3b1SBarry Smith }
122053acd3b1SBarry Smith 
122153acd3b1SBarry Smith #undef __FUNCT__
1222e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsBoolGroup_Private"
122353acd3b1SBarry Smith /*@C
1224acfcf0e5SJed Brown      PetscOptionsBoolGroup - One in a series of logical queries on the options database for
1225d5649816SBarry Smith        which at most a single value can be true.
122653acd3b1SBarry Smith 
12273f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
122853acd3b1SBarry Smith 
122953acd3b1SBarry Smith    Input Parameters:
123053acd3b1SBarry Smith +  opt - option name
123153acd3b1SBarry Smith .  text - short string that describes the option
123253acd3b1SBarry Smith -  man - manual page with additional information on option
123353acd3b1SBarry Smith 
123453acd3b1SBarry Smith    Output Parameter:
123553acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
123653acd3b1SBarry Smith 
123753acd3b1SBarry Smith    Level: intermediate
123853acd3b1SBarry Smith 
123953acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
124053acd3b1SBarry Smith 
1241acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin() and preceded a PetscOptionsBoolGroupEnd()
124253acd3b1SBarry Smith 
124353acd3b1SBarry Smith     Concepts: options database^logical group
124453acd3b1SBarry Smith 
1245c5929fdfSBarry Smith .seealso: PetscOptionsGetInt(NULL,), PetscOptionsGetReal(),
1246acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
124753acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
124853acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1249acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1250a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
125153acd3b1SBarry Smith @*/
12524416b707SBarry Smith PetscErrorCode  PetscOptionsBoolGroup_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool  *flg)
125353acd3b1SBarry Smith {
125453acd3b1SBarry Smith   PetscErrorCode  ierr;
12554416b707SBarry Smith   PetscOptionItem amsopt;
125653acd3b1SBarry Smith 
125753acd3b1SBarry Smith   PetscFunctionBegin;
1258e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
12594416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1260ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1261a297a907SKarl Rupp 
1262ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
12631ae3d29cSBarry Smith   }
126417326d04SJed Brown   *flg = PETSC_FALSE;
1265c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,flg,NULL);CHKERRQ(ierr);
1266e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1267e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"    -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
126853acd3b1SBarry Smith   }
126953acd3b1SBarry Smith   PetscFunctionReturn(0);
127053acd3b1SBarry Smith }
127153acd3b1SBarry Smith 
127253acd3b1SBarry Smith #undef __FUNCT__
1273e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsBoolGroupEnd_Private"
127453acd3b1SBarry Smith /*@C
1275acfcf0e5SJed Brown      PetscOptionsBoolGroupEnd - Last in a series of logical queries on the options database for
1276d5649816SBarry Smith        which at most a single value can be true.
127753acd3b1SBarry Smith 
12783f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
127953acd3b1SBarry Smith 
128053acd3b1SBarry Smith    Input Parameters:
128153acd3b1SBarry Smith +  opt - option name
128253acd3b1SBarry Smith .  text - short string that describes the option
128353acd3b1SBarry Smith -  man - manual page with additional information on option
128453acd3b1SBarry Smith 
128553acd3b1SBarry Smith    Output Parameter:
128653acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
128753acd3b1SBarry Smith 
128853acd3b1SBarry Smith    Level: intermediate
128953acd3b1SBarry Smith 
129053acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
129153acd3b1SBarry Smith 
1292acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin()
129353acd3b1SBarry Smith 
129453acd3b1SBarry Smith     Concepts: options database^logical group
129553acd3b1SBarry Smith 
1296c5929fdfSBarry Smith .seealso: PetscOptionsGetInt(NULL,), PetscOptionsGetReal(),
1297acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
129853acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
129953acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1300acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1301a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
130253acd3b1SBarry Smith @*/
13034416b707SBarry Smith PetscErrorCode  PetscOptionsBoolGroupEnd_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool  *flg)
130453acd3b1SBarry Smith {
130553acd3b1SBarry Smith   PetscErrorCode  ierr;
13064416b707SBarry Smith   PetscOptionItem amsopt;
130753acd3b1SBarry Smith 
130853acd3b1SBarry Smith   PetscFunctionBegin;
1309e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
13104416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1311ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1312a297a907SKarl Rupp 
1313ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
13141ae3d29cSBarry Smith   }
131517326d04SJed Brown   *flg = PETSC_FALSE;
1316c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,flg,NULL);CHKERRQ(ierr);
1317e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1318e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"    -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
131953acd3b1SBarry Smith   }
132053acd3b1SBarry Smith   PetscFunctionReturn(0);
132153acd3b1SBarry Smith }
132253acd3b1SBarry Smith 
132353acd3b1SBarry Smith #undef __FUNCT__
1324e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsBool_Private"
132553acd3b1SBarry Smith /*@C
1326acfcf0e5SJed Brown    PetscOptionsBool - Determines if a particular option is in the database with a true or false
132753acd3b1SBarry Smith 
13283f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
132953acd3b1SBarry Smith 
133053acd3b1SBarry Smith    Input Parameters:
133153acd3b1SBarry Smith +  opt - option name
133253acd3b1SBarry Smith .  text - short string that describes the option
1333868c398cSBarry Smith .  man - manual page with additional information on option
133494ae4db5SBarry Smith -  currentvalue - the current value
133553acd3b1SBarry Smith 
133653acd3b1SBarry Smith    Output Parameter:
133753acd3b1SBarry Smith .  flg - PETSC_TRUE or PETSC_FALSE
133853acd3b1SBarry Smith .  set - PETSC_TRUE if found, else PETSC_FALSE
133953acd3b1SBarry Smith 
134053acd3b1SBarry Smith    Level: beginner
134153acd3b1SBarry Smith 
134253acd3b1SBarry Smith    Concepts: options database^logical
134353acd3b1SBarry Smith 
134453acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
134553acd3b1SBarry Smith 
1346c5929fdfSBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(NULL,),
1347acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
1348acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
134953acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
135053acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1351acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1352a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
135353acd3b1SBarry Smith @*/
13544416b707SBarry Smith PetscErrorCode  PetscOptionsBool_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool currentvalue,PetscBool  *flg,PetscBool  *set)
135553acd3b1SBarry Smith {
135653acd3b1SBarry Smith   PetscErrorCode  ierr;
1357ace3abfcSBarry Smith   PetscBool       iset;
13584416b707SBarry Smith   PetscOptionItem amsopt;
135953acd3b1SBarry Smith 
136053acd3b1SBarry Smith   PetscFunctionBegin;
1361e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
13624416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1363ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1364a297a907SKarl Rupp 
136594ae4db5SBarry Smith     *(PetscBool*)amsopt->data = currentvalue;
1366af6d86caSBarry Smith   }
1367c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,flg,&iset);CHKERRQ(ierr);
136853acd3b1SBarry Smith   if (set) *set = iset;
13691a1499c8SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
137094ae4db5SBarry Smith     const char *v = PetscBools[currentvalue];
13711a1499c8SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s: <%s> %s (%s)\n",PetscOptionsObject->prefix?PetscOptionsObject->prefix:"",opt+1,v,text,ManSection(man));CHKERRQ(ierr);
137253acd3b1SBarry Smith   }
137353acd3b1SBarry Smith   PetscFunctionReturn(0);
137453acd3b1SBarry Smith }
137553acd3b1SBarry Smith 
137653acd3b1SBarry Smith #undef __FUNCT__
1377e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsRealArray_Private"
137853acd3b1SBarry Smith /*@C
137953acd3b1SBarry Smith    PetscOptionsRealArray - Gets an array of double values for a particular
138053acd3b1SBarry Smith    option in the database. The values must be separated with commas with
138153acd3b1SBarry Smith    no intervening spaces.
138253acd3b1SBarry Smith 
13833f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
138453acd3b1SBarry Smith 
138553acd3b1SBarry Smith    Input Parameters:
138653acd3b1SBarry Smith +  opt - the option one is seeking
138753acd3b1SBarry Smith .  text - short string describing option
138853acd3b1SBarry Smith .  man - manual page for option
138953acd3b1SBarry Smith -  nmax - maximum number of values
139053acd3b1SBarry Smith 
139153acd3b1SBarry Smith    Output Parameter:
139253acd3b1SBarry Smith +  value - location to copy values
139353acd3b1SBarry Smith .  nmax - actual number of values found
139453acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
139553acd3b1SBarry Smith 
139653acd3b1SBarry Smith    Level: beginner
139753acd3b1SBarry Smith 
139853acd3b1SBarry Smith    Notes:
139953acd3b1SBarry Smith    The user should pass in an array of doubles
140053acd3b1SBarry Smith 
140153acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
140253acd3b1SBarry Smith 
140353acd3b1SBarry Smith    Concepts: options database^array of strings
140453acd3b1SBarry Smith 
1405c5929fdfSBarry Smith .seealso: PetscOptionsGetInt(NULL,), PetscOptionsGetReal(),
1406acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
140753acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
140853acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1409acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1410a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
141153acd3b1SBarry Smith @*/
14124416b707SBarry Smith PetscErrorCode PetscOptionsRealArray_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscReal value[],PetscInt *n,PetscBool  *set)
141353acd3b1SBarry Smith {
141453acd3b1SBarry Smith   PetscErrorCode  ierr;
141553acd3b1SBarry Smith   PetscInt        i;
14164416b707SBarry Smith   PetscOptionItem amsopt;
141753acd3b1SBarry Smith 
141853acd3b1SBarry Smith   PetscFunctionBegin;
1419e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
1420e26ddf31SBarry Smith     PetscReal *vals;
1421e26ddf31SBarry Smith 
14224416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_REAL_ARRAY,&amsopt);CHKERRQ(ierr);
1423e55864a3SBarry Smith     ierr = PetscMalloc((*n)*sizeof(PetscReal),&amsopt->data);CHKERRQ(ierr);
1424e26ddf31SBarry Smith     vals = (PetscReal*)amsopt->data;
1425e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1426e26ddf31SBarry Smith     amsopt->arraylength = *n;
1427e26ddf31SBarry Smith   }
1428c5929fdfSBarry Smith   ierr = PetscOptionsGetRealArray(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,n,set);CHKERRQ(ierr);
1429e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1430a519f713SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%g",PetscOptionsObject->prefix?PetscOptionsObject->prefix:"",opt+1,(double)value[0]);CHKERRQ(ierr);
143153acd3b1SBarry Smith     for (i=1; i<*n; i++) {
1432e55864a3SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%g",(double)value[i]);CHKERRQ(ierr);
143353acd3b1SBarry Smith     }
1434e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
143553acd3b1SBarry Smith   }
143653acd3b1SBarry Smith   PetscFunctionReturn(0);
143753acd3b1SBarry Smith }
143853acd3b1SBarry Smith 
1439050cccc3SHong Zhang #undef __FUNCT__
1440050cccc3SHong Zhang #define __FUNCT__ "PetscOptionsScalarArray_Private"
1441050cccc3SHong Zhang /*@C
1442050cccc3SHong Zhang    PetscOptionsScalarArray - Gets an array of Scalar values for a particular
1443050cccc3SHong Zhang    option in the database. The values must be separated with commas with
1444050cccc3SHong Zhang    no intervening spaces.
1445050cccc3SHong Zhang 
1446050cccc3SHong Zhang    Logically Collective on the communicator passed in PetscOptionsBegin()
1447050cccc3SHong Zhang 
1448050cccc3SHong Zhang    Input Parameters:
1449050cccc3SHong Zhang +  opt - the option one is seeking
1450050cccc3SHong Zhang .  text - short string describing option
1451050cccc3SHong Zhang .  man - manual page for option
1452050cccc3SHong Zhang -  nmax - maximum number of values
1453050cccc3SHong Zhang 
1454050cccc3SHong Zhang    Output Parameter:
1455050cccc3SHong Zhang +  value - location to copy values
1456050cccc3SHong Zhang .  nmax - actual number of values found
1457050cccc3SHong Zhang -  set - PETSC_TRUE if found, else PETSC_FALSE
1458050cccc3SHong Zhang 
1459050cccc3SHong Zhang    Level: beginner
1460050cccc3SHong Zhang 
1461050cccc3SHong Zhang    Notes:
1462050cccc3SHong Zhang    The user should pass in an array of doubles
1463050cccc3SHong Zhang 
1464050cccc3SHong Zhang    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
1465050cccc3SHong Zhang 
1466050cccc3SHong Zhang    Concepts: options database^array of strings
1467050cccc3SHong Zhang 
1468c5929fdfSBarry Smith .seealso: PetscOptionsGetInt(NULL,), PetscOptionsGetReal(),
1469050cccc3SHong Zhang           PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
1470050cccc3SHong Zhang           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1471050cccc3SHong Zhang           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1472050cccc3SHong Zhang           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1473050cccc3SHong Zhang           PetscOptionsFList(), PetscOptionsEList()
1474050cccc3SHong Zhang @*/
14754416b707SBarry Smith PetscErrorCode PetscOptionsScalarArray_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscScalar value[],PetscInt *n,PetscBool  *set)
1476050cccc3SHong Zhang {
1477050cccc3SHong Zhang   PetscErrorCode  ierr;
1478050cccc3SHong Zhang   PetscInt        i;
14794416b707SBarry Smith   PetscOptionItem amsopt;
1480050cccc3SHong Zhang 
1481050cccc3SHong Zhang   PetscFunctionBegin;
1482050cccc3SHong Zhang   if (!PetscOptionsObject->count) {
1483050cccc3SHong Zhang     PetscScalar *vals;
1484050cccc3SHong Zhang 
14854416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_SCALAR_ARRAY,&amsopt);CHKERRQ(ierr);
1486050cccc3SHong Zhang     ierr = PetscMalloc((*n)*sizeof(PetscScalar),&amsopt->data);CHKERRQ(ierr);
1487050cccc3SHong Zhang     vals = (PetscScalar*)amsopt->data;
1488050cccc3SHong Zhang     for (i=0; i<*n; i++) vals[i] = value[i];
1489050cccc3SHong Zhang     amsopt->arraylength = *n;
1490050cccc3SHong Zhang   }
1491c5929fdfSBarry Smith   ierr = PetscOptionsGetScalarArray(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,n,set);CHKERRQ(ierr);
1492050cccc3SHong Zhang   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
149395f3a755SJose 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);
1494050cccc3SHong Zhang     for (i=1; i<*n; i++) {
149595f3a755SJose E. Roman       ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%g+%gi",(double)PetscRealPart(value[i]),(double)PetscImaginaryPart(value[i]));CHKERRQ(ierr);
1496050cccc3SHong Zhang     }
1497050cccc3SHong Zhang     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
1498050cccc3SHong Zhang   }
1499050cccc3SHong Zhang   PetscFunctionReturn(0);
1500050cccc3SHong Zhang }
150153acd3b1SBarry Smith 
150253acd3b1SBarry Smith #undef __FUNCT__
1503e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsIntArray_Private"
150453acd3b1SBarry Smith /*@C
150553acd3b1SBarry Smith    PetscOptionsIntArray - Gets an array of integers for a particular
1506b32a342fSShri Abhyankar    option in the database.
150753acd3b1SBarry Smith 
15083f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
150953acd3b1SBarry Smith 
151053acd3b1SBarry Smith    Input Parameters:
151153acd3b1SBarry Smith +  opt - the option one is seeking
151253acd3b1SBarry Smith .  text - short string describing option
151353acd3b1SBarry Smith .  man - manual page for option
1514f8a50e2bSBarry Smith -  n - maximum number of values
151553acd3b1SBarry Smith 
151653acd3b1SBarry Smith    Output Parameter:
151753acd3b1SBarry Smith +  value - location to copy values
1518f8a50e2bSBarry Smith .  n - actual number of values found
151953acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
152053acd3b1SBarry Smith 
152153acd3b1SBarry Smith    Level: beginner
152253acd3b1SBarry Smith 
152353acd3b1SBarry Smith    Notes:
1524b32a342fSShri Abhyankar    The array can be passed as
1525bebe2cf6SSatish Balay    a comma separated list:                                 0,1,2,3,4,5,6,7
15260fd488f5SShri Abhyankar    a range (start-end+1):                                  0-8
15270fd488f5SShri Abhyankar    a range with given increment (start-end+1:inc):         0-7:2
1528bebe2cf6SSatish Balay    a combination of values and ranges separated by commas: 0,1-8,8-15:2
1529b32a342fSShri Abhyankar 
1530b32a342fSShri Abhyankar    There must be no intervening spaces between the values.
153153acd3b1SBarry Smith 
153253acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
153353acd3b1SBarry Smith 
1534b32a342fSShri Abhyankar    Concepts: options database^array of ints
153553acd3b1SBarry Smith 
1536c5929fdfSBarry Smith .seealso: PetscOptionsGetInt(NULL,), PetscOptionsGetReal(),
1537acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
153853acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
153953acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1540acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1541a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsRealArray()
154253acd3b1SBarry Smith @*/
15434416b707SBarry Smith PetscErrorCode  PetscOptionsIntArray_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscInt value[],PetscInt *n,PetscBool  *set)
154453acd3b1SBarry Smith {
154553acd3b1SBarry Smith   PetscErrorCode ierr;
154653acd3b1SBarry Smith   PetscInt        i;
15474416b707SBarry Smith   PetscOptionItem amsopt;
154853acd3b1SBarry Smith 
154953acd3b1SBarry Smith   PetscFunctionBegin;
1550e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
1551e26ddf31SBarry Smith     PetscInt *vals;
1552e26ddf31SBarry Smith 
15534416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_INT_ARRAY,&amsopt);CHKERRQ(ierr);
1554854ce69bSBarry Smith     ierr = PetscMalloc1(*n,(PetscInt**)&amsopt->data);CHKERRQ(ierr);
1555e26ddf31SBarry Smith     vals = (PetscInt*)amsopt->data;
1556e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1557e26ddf31SBarry Smith     amsopt->arraylength = *n;
1558e26ddf31SBarry Smith   }
1559c5929fdfSBarry Smith   ierr = PetscOptionsGetIntArray(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,n,set);CHKERRQ(ierr);
1560e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1561e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%d",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,value[0]);CHKERRQ(ierr);
156253acd3b1SBarry Smith     for (i=1; i<*n; i++) {
1563e55864a3SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%d",value[i]);CHKERRQ(ierr);
156453acd3b1SBarry Smith     }
1565e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
156653acd3b1SBarry Smith   }
156753acd3b1SBarry Smith   PetscFunctionReturn(0);
156853acd3b1SBarry Smith }
156953acd3b1SBarry Smith 
157053acd3b1SBarry Smith #undef __FUNCT__
1571e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsStringArray_Private"
157253acd3b1SBarry Smith /*@C
157353acd3b1SBarry Smith    PetscOptionsStringArray - Gets an array of string values for a particular
157453acd3b1SBarry Smith    option in the database. The values must be separated with commas with
157553acd3b1SBarry Smith    no intervening spaces.
157653acd3b1SBarry Smith 
15773f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
157853acd3b1SBarry Smith 
157953acd3b1SBarry Smith    Input Parameters:
158053acd3b1SBarry Smith +  opt - the option one is seeking
158153acd3b1SBarry Smith .  text - short string describing option
158253acd3b1SBarry Smith .  man - manual page for option
158353acd3b1SBarry Smith -  nmax - maximum number of strings
158453acd3b1SBarry Smith 
158553acd3b1SBarry Smith    Output Parameter:
158653acd3b1SBarry Smith +  value - location to copy strings
158753acd3b1SBarry Smith .  nmax - actual number of strings found
158853acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
158953acd3b1SBarry Smith 
159053acd3b1SBarry Smith    Level: beginner
159153acd3b1SBarry Smith 
159253acd3b1SBarry Smith    Notes:
159353acd3b1SBarry Smith    The user should pass in an array of pointers to char, to hold all the
159453acd3b1SBarry Smith    strings returned by this function.
159553acd3b1SBarry Smith 
159653acd3b1SBarry Smith    The user is responsible for deallocating the strings that are
159753acd3b1SBarry Smith    returned. The Fortran interface for this routine is not supported.
159853acd3b1SBarry Smith 
159953acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
160053acd3b1SBarry Smith 
160153acd3b1SBarry Smith    Concepts: options database^array of strings
160253acd3b1SBarry Smith 
1603c5929fdfSBarry Smith .seealso: PetscOptionsGetInt(NULL,), PetscOptionsGetReal(),
1604acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
160553acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
160653acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1607acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1608a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
160953acd3b1SBarry Smith @*/
16104416b707SBarry Smith PetscErrorCode  PetscOptionsStringArray_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],char *value[],PetscInt *nmax,PetscBool  *set)
161153acd3b1SBarry Smith {
161253acd3b1SBarry Smith   PetscErrorCode  ierr;
16134416b707SBarry Smith   PetscOptionItem amsopt;
161453acd3b1SBarry Smith 
161553acd3b1SBarry Smith   PetscFunctionBegin;
1616e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
16174416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_STRING_ARRAY,&amsopt);CHKERRQ(ierr);
1618854ce69bSBarry Smith     ierr = PetscMalloc1(*nmax,(char**)&amsopt->data);CHKERRQ(ierr);
1619a297a907SKarl Rupp 
16201ae3d29cSBarry Smith     amsopt->arraylength = *nmax;
16211ae3d29cSBarry Smith   }
1622c5929fdfSBarry Smith   ierr = PetscOptionsGetStringArray(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,nmax,set);CHKERRQ(ierr);
1623e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1624e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <string1,string2,...>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
162553acd3b1SBarry Smith   }
162653acd3b1SBarry Smith   PetscFunctionReturn(0);
162753acd3b1SBarry Smith }
162853acd3b1SBarry Smith 
1629e2446a98SMatthew Knepley #undef __FUNCT__
1630e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsBoolArray_Private"
1631e2446a98SMatthew Knepley /*@C
1632acfcf0e5SJed Brown    PetscOptionsBoolArray - Gets an array of logical values (true or false) for a particular
1633e2446a98SMatthew Knepley    option in the database. The values must be separated with commas with
1634e2446a98SMatthew Knepley    no intervening spaces.
1635e2446a98SMatthew Knepley 
16363f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
1637e2446a98SMatthew Knepley 
1638e2446a98SMatthew Knepley    Input Parameters:
1639e2446a98SMatthew Knepley +  opt - the option one is seeking
1640e2446a98SMatthew Knepley .  text - short string describing option
1641e2446a98SMatthew Knepley .  man - manual page for option
1642e2446a98SMatthew Knepley -  nmax - maximum number of values
1643e2446a98SMatthew Knepley 
1644e2446a98SMatthew Knepley    Output Parameter:
1645e2446a98SMatthew Knepley +  value - location to copy values
1646e2446a98SMatthew Knepley .  nmax - actual number of values found
1647e2446a98SMatthew Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
1648e2446a98SMatthew Knepley 
1649e2446a98SMatthew Knepley    Level: beginner
1650e2446a98SMatthew Knepley 
1651e2446a98SMatthew Knepley    Notes:
1652e2446a98SMatthew Knepley    The user should pass in an array of doubles
1653e2446a98SMatthew Knepley 
1654e2446a98SMatthew Knepley    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
1655e2446a98SMatthew Knepley 
1656e2446a98SMatthew Knepley    Concepts: options database^array of strings
1657e2446a98SMatthew Knepley 
1658c5929fdfSBarry Smith .seealso: PetscOptionsGetInt(NULL,), PetscOptionsGetReal(),
1659acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
1660e2446a98SMatthew Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1661e2446a98SMatthew Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1662acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1663a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
1664e2446a98SMatthew Knepley @*/
16654416b707SBarry Smith PetscErrorCode  PetscOptionsBoolArray_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool value[],PetscInt *n,PetscBool *set)
1666e2446a98SMatthew Knepley {
1667e2446a98SMatthew Knepley   PetscErrorCode   ierr;
1668e2446a98SMatthew Knepley   PetscInt         i;
16694416b707SBarry Smith   PetscOptionItem  amsopt;
1670e2446a98SMatthew Knepley 
1671e2446a98SMatthew Knepley   PetscFunctionBegin;
1672e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
1673ace3abfcSBarry Smith     PetscBool *vals;
16741ae3d29cSBarry Smith 
16754416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL_ARRAY,&amsopt);CHKERRQ(ierr);
16761a1499c8SBarry Smith     ierr = PetscMalloc1(*n,(PetscBool**)&amsopt->data);CHKERRQ(ierr);
1677ace3abfcSBarry Smith     vals = (PetscBool*)amsopt->data;
16781ae3d29cSBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
16791ae3d29cSBarry Smith     amsopt->arraylength = *n;
16801ae3d29cSBarry Smith   }
1681c5929fdfSBarry Smith   ierr = PetscOptionsGetBoolArray(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,n,set);CHKERRQ(ierr);
1682e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1683e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%d",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,value[0]);CHKERRQ(ierr);
1684e2446a98SMatthew Knepley     for (i=1; i<*n; i++) {
1685e55864a3SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%d",value[i]);CHKERRQ(ierr);
1686e2446a98SMatthew Knepley     }
1687e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
1688e2446a98SMatthew Knepley   }
1689e2446a98SMatthew Knepley   PetscFunctionReturn(0);
1690e2446a98SMatthew Knepley }
1691e2446a98SMatthew Knepley 
16928cc676e6SMatthew G Knepley #undef __FUNCT__
1693e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsViewer_Private"
16948cc676e6SMatthew G Knepley /*@C
1695d1da0b69SBarry Smith    PetscOptionsViewer - Gets a viewer appropriate for the type indicated by the user
16968cc676e6SMatthew G Knepley 
16978cc676e6SMatthew G Knepley    Logically Collective on the communicator passed in PetscOptionsBegin()
16988cc676e6SMatthew G Knepley 
16998cc676e6SMatthew G Knepley    Input Parameters:
17008cc676e6SMatthew G Knepley +  opt - option name
17018cc676e6SMatthew G Knepley .  text - short string that describes the option
17028cc676e6SMatthew G Knepley -  man - manual page with additional information on option
17038cc676e6SMatthew G Knepley 
17048cc676e6SMatthew G Knepley    Output Parameter:
17058cc676e6SMatthew G Knepley +  viewer - the viewer
17068cc676e6SMatthew G Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
17078cc676e6SMatthew G Knepley 
17088cc676e6SMatthew G Knepley    Level: beginner
17098cc676e6SMatthew G Knepley 
17108cc676e6SMatthew G Knepley    Concepts: options database^has int
17118cc676e6SMatthew G Knepley 
17128cc676e6SMatthew G Knepley    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
17138cc676e6SMatthew G Knepley 
17145a7113b9SPatrick Sanan    See PetscOptionsGetViewer() for the format of the supplied viewer and its options
17158cc676e6SMatthew G Knepley 
1716c5929fdfSBarry Smith .seealso: PetscOptionsGetViewer(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(NULL,),
17178cc676e6SMatthew G Knepley           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
17188cc676e6SMatthew G Knepley           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
17198cc676e6SMatthew G Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
17208cc676e6SMatthew G Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
17218cc676e6SMatthew G Knepley           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1722a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
17238cc676e6SMatthew G Knepley @*/
17244416b707SBarry Smith PetscErrorCode  PetscOptionsViewer_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscViewer *viewer,PetscViewerFormat *format,PetscBool  *set)
17258cc676e6SMatthew G Knepley {
17268cc676e6SMatthew G Knepley   PetscErrorCode  ierr;
17274416b707SBarry Smith   PetscOptionItem amsopt;
17288cc676e6SMatthew G Knepley 
17298cc676e6SMatthew G Knepley   PetscFunctionBegin;
17301a1499c8SBarry Smith   if (!PetscOptionsObject->count) {
17314416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
173264facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
17335b02f95dSBarry Smith     ierr = PetscStrdup("",(char**)&amsopt->data);CHKERRQ(ierr);
17348cc676e6SMatthew G Knepley   }
1735e55864a3SBarry Smith   ierr = PetscOptionsGetViewer(PetscOptionsObject->comm,PetscOptionsObject->prefix,opt,viewer,format,set);CHKERRQ(ierr);
1736e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1737e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%s>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,"",text,ManSection(man));CHKERRQ(ierr);
17388cc676e6SMatthew G Knepley   }
17398cc676e6SMatthew G Knepley   PetscFunctionReturn(0);
17408cc676e6SMatthew G Knepley }
17418cc676e6SMatthew G Knepley 
174253acd3b1SBarry Smith 
174353acd3b1SBarry Smith #undef __FUNCT__
174453acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsHead"
174553acd3b1SBarry Smith /*@C
1746b52f573bSBarry Smith      PetscOptionsHead - Puts a heading before listing any more published options. Used, for example,
174753acd3b1SBarry Smith             in KSPSetFromOptions_GMRES().
174853acd3b1SBarry Smith 
17493f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
175053acd3b1SBarry Smith 
175153acd3b1SBarry Smith    Input Parameter:
175253acd3b1SBarry Smith .   head - the heading text
175353acd3b1SBarry Smith 
175453acd3b1SBarry Smith 
175553acd3b1SBarry Smith    Level: intermediate
175653acd3b1SBarry Smith 
175753acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
175853acd3b1SBarry Smith 
1759b52f573bSBarry Smith           Can be followed by a call to PetscOptionsTail() in the same function.
176053acd3b1SBarry Smith 
176153acd3b1SBarry Smith    Concepts: options database^subheading
176253acd3b1SBarry Smith 
1763c5929fdfSBarry Smith .seealso: PetscOptionsGetInt(NULL,), PetscOptionsGetReal(),
1764acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
176553acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
176653acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1767acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1768a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
176953acd3b1SBarry Smith @*/
17704416b707SBarry Smith PetscErrorCode  PetscOptionsHead(PetscOptionItems *PetscOptionsObject,const char head[])
177153acd3b1SBarry Smith {
177253acd3b1SBarry Smith   PetscErrorCode ierr;
177353acd3b1SBarry Smith 
177453acd3b1SBarry Smith   PetscFunctionBegin;
1775e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1776e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  %s\n",head);CHKERRQ(ierr);
177753acd3b1SBarry Smith   }
177853acd3b1SBarry Smith   PetscFunctionReturn(0);
177953acd3b1SBarry Smith }
178053acd3b1SBarry Smith 
178153acd3b1SBarry Smith 
178253acd3b1SBarry Smith 
178353acd3b1SBarry Smith 
178453acd3b1SBarry Smith 
178553acd3b1SBarry Smith 
1786