xref: /petsc/src/sys/objects/aoptions.c (revision 989712b93983239785d01a070533fb7e89598ec3)
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 /*
2353acd3b1SBarry Smith     Handles setting up the data structure in a call to PetscOptionsBegin()
2453acd3b1SBarry Smith */
254416b707SBarry Smith PetscErrorCode PetscOptionsBegin_Private(PetscOptionItems *PetscOptionsObject,MPI_Comm comm,const char prefix[],const char title[],const char mansec[])
2653acd3b1SBarry Smith {
2753acd3b1SBarry Smith   PetscErrorCode ierr;
2853acd3b1SBarry Smith 
2953acd3b1SBarry Smith   PetscFunctionBegin;
300eb63584SBarry Smith   if (!PetscOptionsObject->alreadyprinted) {
319de0f6ecSBarry Smith     if (!PetscOptionsHelpPrintedSingleton) {
329de0f6ecSBarry Smith       ierr = PetscOptionsHelpPrintedCreate(&PetscOptionsHelpPrintedSingleton);CHKERRQ(ierr);
339de0f6ecSBarry Smith     }
349de0f6ecSBarry Smith     ierr = PetscOptionsHelpPrintedCheck(PetscOptionsHelpPrintedSingleton,prefix,title,&PetscOptionsObject->alreadyprinted);CHKERRQ(ierr);
350eb63584SBarry Smith   }
36e55864a3SBarry Smith   PetscOptionsObject->next          = 0;
37e55864a3SBarry Smith   PetscOptionsObject->comm          = comm;
38e55864a3SBarry Smith   PetscOptionsObject->changedmethod = PETSC_FALSE;
39a297a907SKarl Rupp 
40e55864a3SBarry Smith   ierr = PetscStrallocpy(prefix,&PetscOptionsObject->prefix);CHKERRQ(ierr);
41e55864a3SBarry Smith   ierr = PetscStrallocpy(title,&PetscOptionsObject->title);CHKERRQ(ierr);
4253acd3b1SBarry Smith 
43c5929fdfSBarry Smith   ierr = PetscOptionsHasName(PetscOptionsObject->options,NULL,"-help",&PetscOptionsObject->printhelp);CHKERRQ(ierr);
44e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1) {
45e55864a3SBarry Smith     if (!PetscOptionsObject->alreadyprinted) {
4653acd3b1SBarry Smith       ierr = (*PetscHelpPrintf)(comm,"%s -------------------------------------------------\n",title);CHKERRQ(ierr);
4753acd3b1SBarry Smith     }
4861b37b28SSatish Balay   }
4953acd3b1SBarry Smith   PetscFunctionReturn(0);
5053acd3b1SBarry Smith }
5153acd3b1SBarry Smith 
523194b578SJed Brown /*
533194b578SJed Brown     Handles setting up the data structure in a call to PetscObjectOptionsBegin()
543194b578SJed Brown */
554416b707SBarry Smith PetscErrorCode PetscObjectOptionsBegin_Private(PetscOptionItems *PetscOptionsObject,PetscObject obj)
563194b578SJed Brown {
573194b578SJed Brown   PetscErrorCode ierr;
583194b578SJed Brown   char           title[256];
593194b578SJed Brown   PetscBool      flg;
603194b578SJed Brown 
613194b578SJed Brown   PetscFunctionBegin;
623194b578SJed Brown   PetscValidHeader(obj,1);
63e55864a3SBarry Smith   PetscOptionsObject->object         = obj;
64e55864a3SBarry Smith   PetscOptionsObject->alreadyprinted = obj->optionsprinted;
65a297a907SKarl Rupp 
663194b578SJed Brown   ierr = PetscStrcmp(obj->description,obj->class_name,&flg);CHKERRQ(ierr);
673194b578SJed Brown   if (flg) {
688caf3d72SBarry Smith     ierr = PetscSNPrintf(title,sizeof(title),"%s options",obj->class_name);CHKERRQ(ierr);
693194b578SJed Brown   } else {
708caf3d72SBarry Smith     ierr = PetscSNPrintf(title,sizeof(title),"%s (%s) options",obj->description,obj->class_name);CHKERRQ(ierr);
713194b578SJed Brown   }
72e55864a3SBarry Smith   ierr = PetscOptionsBegin_Private(PetscOptionsObject,obj->comm,obj->prefix,title,obj->mansec);CHKERRQ(ierr);
733194b578SJed Brown   PetscFunctionReturn(0);
743194b578SJed Brown }
753194b578SJed Brown 
7653acd3b1SBarry Smith /*
7753acd3b1SBarry Smith      Handles adding another option to the list of options within this particular PetscOptionsBegin() PetscOptionsEnd()
7853acd3b1SBarry Smith */
794416b707SBarry Smith static int PetscOptionItemCreate_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscOptionType t,PetscOptionItem *amsopt)
8053acd3b1SBarry Smith {
8153acd3b1SBarry Smith   int             ierr;
824416b707SBarry Smith   PetscOptionItem next;
833be6e4c3SJed Brown   PetscBool       valid;
8453acd3b1SBarry Smith 
8553acd3b1SBarry Smith   PetscFunctionBegin;
863be6e4c3SJed Brown   ierr = PetscOptionsValidKey(opt,&valid);CHKERRQ(ierr);
873be6e4c3SJed Brown   if (!valid) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"The option '%s' is not a valid key",opt);
883be6e4c3SJed Brown 
89b00a9115SJed Brown   ierr            = PetscNew(amsopt);CHKERRQ(ierr);
9053acd3b1SBarry Smith   (*amsopt)->next = 0;
9153acd3b1SBarry Smith   (*amsopt)->set  = PETSC_FALSE;
926356e834SBarry Smith   (*amsopt)->type = t;
9353acd3b1SBarry Smith   (*amsopt)->data = 0;
9461b37b28SSatish Balay 
9553acd3b1SBarry Smith   ierr = PetscStrallocpy(text,&(*amsopt)->text);CHKERRQ(ierr);
9653acd3b1SBarry Smith   ierr = PetscStrallocpy(opt,&(*amsopt)->option);CHKERRQ(ierr);
976356e834SBarry Smith   ierr = PetscStrallocpy(man,&(*amsopt)->man);CHKERRQ(ierr);
9853acd3b1SBarry Smith 
99e55864a3SBarry Smith   if (!PetscOptionsObject->next) PetscOptionsObject->next = *amsopt;
100a297a907SKarl Rupp   else {
101e55864a3SBarry Smith     next = PetscOptionsObject->next;
10253acd3b1SBarry Smith     while (next->next) next = next->next;
10353acd3b1SBarry Smith     next->next = *amsopt;
10453acd3b1SBarry Smith   }
10553acd3b1SBarry Smith   PetscFunctionReturn(0);
10653acd3b1SBarry Smith }
10753acd3b1SBarry Smith 
108aee2cecaSBarry Smith /*
1093fc1eb6aSBarry Smith     PetscScanString -  Gets user input via stdin from process and broadcasts to all processes
1103fc1eb6aSBarry Smith 
1113fc1eb6aSBarry Smith     Collective on MPI_Comm
1123fc1eb6aSBarry Smith 
1133fc1eb6aSBarry Smith    Input Parameters:
1143fc1eb6aSBarry Smith +     commm - communicator for the broadcast, must be PETSC_COMM_WORLD
1153fc1eb6aSBarry Smith .     n - length of the string, must be the same on all processes
1163fc1eb6aSBarry Smith -     str - location to store input
117aee2cecaSBarry Smith 
118aee2cecaSBarry Smith     Bugs:
119aee2cecaSBarry Smith .   Assumes process 0 of the given communicator has access to stdin
120aee2cecaSBarry Smith 
121aee2cecaSBarry Smith */
1223fc1eb6aSBarry Smith static PetscErrorCode PetscScanString(MPI_Comm comm,size_t n,char str[])
123aee2cecaSBarry Smith {
124330cf3c9SBarry Smith   size_t         i;
125aee2cecaSBarry Smith   char           c;
1263fc1eb6aSBarry Smith   PetscMPIInt    rank,nm;
127aee2cecaSBarry Smith   PetscErrorCode ierr;
128aee2cecaSBarry Smith 
129aee2cecaSBarry Smith   PetscFunctionBegin;
130aee2cecaSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
131aee2cecaSBarry Smith   if (!rank) {
132aee2cecaSBarry Smith     c = (char) getchar();
133aee2cecaSBarry Smith     i = 0;
134aee2cecaSBarry Smith     while (c != '\n' && i < n-1) {
135aee2cecaSBarry Smith       str[i++] = c;
136aee2cecaSBarry Smith       c = (char)getchar();
137aee2cecaSBarry Smith     }
138aee2cecaSBarry Smith     str[i] = 0;
139aee2cecaSBarry Smith   }
1404dc2109aSBarry Smith   ierr = PetscMPIIntCast(n,&nm);CHKERRQ(ierr);
1413fc1eb6aSBarry Smith   ierr = MPI_Bcast(str,nm,MPI_CHAR,0,comm);CHKERRQ(ierr);
142aee2cecaSBarry Smith   PetscFunctionReturn(0);
143aee2cecaSBarry Smith }
144aee2cecaSBarry Smith 
1455b02f95dSBarry Smith /*
1465b02f95dSBarry Smith     This is needed because certain strings may be freed by SAWs, hence we cannot use PetscStrallocpy()
1475b02f95dSBarry Smith */
1485b02f95dSBarry Smith static PetscErrorCode  PetscStrdup(const char s[],char *t[])
1495b02f95dSBarry Smith {
1505b02f95dSBarry Smith   PetscErrorCode ierr;
1515b02f95dSBarry Smith   size_t         len;
1525b02f95dSBarry Smith   char           *tmp = 0;
1535b02f95dSBarry Smith 
1545b02f95dSBarry Smith   PetscFunctionBegin;
1555b02f95dSBarry Smith   if (s) {
1565b02f95dSBarry Smith     ierr = PetscStrlen(s,&len);CHKERRQ(ierr);
157f416af30SBarry Smith     tmp = (char*) malloc((len+1)*sizeof(char));
1585b02f95dSBarry Smith     if (!tmp) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"No memory to duplicate string");
1595b02f95dSBarry Smith     ierr = PetscStrcpy(tmp,s);CHKERRQ(ierr);
1605b02f95dSBarry Smith   }
1615b02f95dSBarry Smith   *t = tmp;
1625b02f95dSBarry Smith   PetscFunctionReturn(0);
1635b02f95dSBarry Smith }
1645b02f95dSBarry Smith 
1655b02f95dSBarry Smith 
166aee2cecaSBarry Smith /*
1673cc1e11dSBarry Smith     PetscOptionsGetFromTextInput - Presents all the PETSc Options processed by the program so the user may change them at runtime
168aee2cecaSBarry Smith 
169aee2cecaSBarry Smith     Notes: this isn't really practical, it is just to demonstrate the principle
170aee2cecaSBarry Smith 
1717781c08eSBarry Smith     A carriage return indicates no change from the default; but this like -ksp_monitor <stdout>  the default is actually not stdout the default
1727781c08eSBarry Smith     is to do nothing so to get it to use stdout you need to type stdout. This is kind of bug?
1737781c08eSBarry Smith 
174aee2cecaSBarry Smith     Bugs:
1757781c08eSBarry Smith +    All processes must traverse through the exact same set of option queries due to the call to PetscScanString()
1763cc1e11dSBarry Smith .    Internal strings have arbitrary length and string copies are not checked that they fit into string space
177aee2cecaSBarry Smith -    Only works for PetscInt == int, PetscReal == double etc
178aee2cecaSBarry Smith 
1793cc1e11dSBarry Smith     Developer Notes: Normally the GUI that presents the options the user and retrieves the values would be running in a different
1803cc1e11dSBarry Smith      address space and communicating with the PETSc program
1813cc1e11dSBarry Smith 
182aee2cecaSBarry Smith */
1834416b707SBarry Smith PetscErrorCode PetscOptionsGetFromTextInput(PetscOptionItems *PetscOptionsObject)
1846356e834SBarry Smith {
1856356e834SBarry Smith   PetscErrorCode  ierr;
1864416b707SBarry Smith   PetscOptionItem next = PetscOptionsObject->next;
1876356e834SBarry Smith   char            str[512];
1887781c08eSBarry Smith   PetscBool       bid;
189a4404d99SBarry Smith   PetscReal       ir,*valr;
190330cf3c9SBarry Smith   PetscInt        *vald;
191330cf3c9SBarry Smith   size_t          i;
1926356e834SBarry Smith 
1932a409bb0SBarry Smith   PetscFunctionBegin;
194e55864a3SBarry Smith   ierr = (*PetscPrintf)(PETSC_COMM_WORLD,"%s -------------------------------------------------\n",PetscOptionsObject->title);CHKERRQ(ierr);
1956356e834SBarry Smith   while (next) {
1966356e834SBarry Smith     switch (next->type) {
1976356e834SBarry Smith     case OPTION_HEAD:
1986356e834SBarry Smith       break;
199e26ddf31SBarry Smith     case OPTION_INT_ARRAY:
200e55864a3SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",next->option+1);CHKERRQ(ierr);
201e26ddf31SBarry Smith       vald = (PetscInt*) next->data;
202e26ddf31SBarry Smith       for (i=0; i<next->arraylength; i++) {
203e26ddf31SBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"%d",vald[i]);CHKERRQ(ierr);
204e26ddf31SBarry Smith         if (i < next->arraylength-1) {
205e26ddf31SBarry Smith           ierr = PetscPrintf(PETSC_COMM_WORLD,",");CHKERRQ(ierr);
206e26ddf31SBarry Smith         }
207e26ddf31SBarry Smith       }
208e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,">: %s (%s) ",next->text,next->man);CHKERRQ(ierr);
209e26ddf31SBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
210e26ddf31SBarry Smith       if (str[0]) {
211e26ddf31SBarry Smith         PetscToken token;
212e26ddf31SBarry Smith         PetscInt   n=0,nmax = next->arraylength,*dvalue = (PetscInt*)next->data,start,end;
213e26ddf31SBarry Smith         size_t     len;
214e26ddf31SBarry Smith         char       *value;
215ace3abfcSBarry Smith         PetscBool  foundrange;
216e26ddf31SBarry Smith 
217e26ddf31SBarry Smith         next->set = PETSC_TRUE;
218e26ddf31SBarry Smith         value     = str;
219e26ddf31SBarry Smith         ierr      = PetscTokenCreate(value,',',&token);CHKERRQ(ierr);
220e26ddf31SBarry Smith         ierr      = PetscTokenFind(token,&value);CHKERRQ(ierr);
221e26ddf31SBarry Smith         while (n < nmax) {
222e26ddf31SBarry Smith           if (!value) break;
223e26ddf31SBarry Smith 
224e26ddf31SBarry Smith           /* look for form  d-D where d and D are integers */
225e26ddf31SBarry Smith           foundrange = PETSC_FALSE;
226e26ddf31SBarry Smith           ierr       = PetscStrlen(value,&len);CHKERRQ(ierr);
227e26ddf31SBarry Smith           if (value[0] == '-') i=2;
228e26ddf31SBarry Smith           else i=1;
229330cf3c9SBarry Smith           for (;i<len; i++) {
230e26ddf31SBarry Smith             if (value[i] == '-') {
231e32f2f54SBarry Smith               if (i == len-1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Error in %D-th array entry %s\n",n,value);
232e26ddf31SBarry Smith               value[i] = 0;
233cfbddea1SSatish Balay               ierr     = PetscOptionsStringToInt(value,&start);CHKERRQ(ierr);
234cfbddea1SSatish Balay               ierr     = PetscOptionsStringToInt(value+i+1,&end);CHKERRQ(ierr);
235e32f2f54SBarry 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);
236e32f2f54SBarry 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);
237e26ddf31SBarry Smith               for (; start<end; start++) {
238e26ddf31SBarry Smith                 *dvalue = start; dvalue++;n++;
239e26ddf31SBarry Smith               }
240e26ddf31SBarry Smith               foundrange = PETSC_TRUE;
241e26ddf31SBarry Smith               break;
242e26ddf31SBarry Smith             }
243e26ddf31SBarry Smith           }
244e26ddf31SBarry Smith           if (!foundrange) {
245cfbddea1SSatish Balay             ierr = PetscOptionsStringToInt(value,dvalue);CHKERRQ(ierr);
246e26ddf31SBarry Smith             dvalue++;
247e26ddf31SBarry Smith             n++;
248e26ddf31SBarry Smith           }
249e26ddf31SBarry Smith           ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
250e26ddf31SBarry Smith         }
2518c74ee41SBarry Smith         ierr = PetscTokenDestroy(&token);CHKERRQ(ierr);
252e26ddf31SBarry Smith       }
253e26ddf31SBarry Smith       break;
254e26ddf31SBarry Smith     case OPTION_REAL_ARRAY:
255e55864a3SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",next->option+1);CHKERRQ(ierr);
256e26ddf31SBarry Smith       valr = (PetscReal*) next->data;
257e26ddf31SBarry Smith       for (i=0; i<next->arraylength; i++) {
258e26ddf31SBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"%g",valr[i]);CHKERRQ(ierr);
259e26ddf31SBarry Smith         if (i < next->arraylength-1) {
260e26ddf31SBarry Smith           ierr = PetscPrintf(PETSC_COMM_WORLD,",");CHKERRQ(ierr);
261e26ddf31SBarry Smith         }
262e26ddf31SBarry Smith       }
263e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,">: %s (%s) ",next->text,next->man);CHKERRQ(ierr);
264e26ddf31SBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
265e26ddf31SBarry Smith       if (str[0]) {
266e26ddf31SBarry Smith         PetscToken token;
267e26ddf31SBarry Smith         PetscInt   n = 0,nmax = next->arraylength;
268e26ddf31SBarry Smith         PetscReal  *dvalue = (PetscReal*)next->data;
269e26ddf31SBarry Smith         char       *value;
270e26ddf31SBarry Smith 
271e26ddf31SBarry Smith         next->set = PETSC_TRUE;
272e26ddf31SBarry Smith         value     = str;
273e26ddf31SBarry Smith         ierr      = PetscTokenCreate(value,',',&token);CHKERRQ(ierr);
274e26ddf31SBarry Smith         ierr      = PetscTokenFind(token,&value);CHKERRQ(ierr);
275e26ddf31SBarry Smith         while (n < nmax) {
276e26ddf31SBarry Smith           if (!value) break;
277cfbddea1SSatish Balay           ierr = PetscOptionsStringToReal(value,dvalue);CHKERRQ(ierr);
278e26ddf31SBarry Smith           dvalue++;
279e26ddf31SBarry Smith           n++;
280e26ddf31SBarry Smith           ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
281e26ddf31SBarry Smith         }
2828c74ee41SBarry Smith         ierr = PetscTokenDestroy(&token);CHKERRQ(ierr);
283e26ddf31SBarry Smith       }
284e26ddf31SBarry Smith       break;
2856356e834SBarry Smith     case OPTION_INT:
286e55864a3SBarry 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);
2873fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
2883fc1eb6aSBarry Smith       if (str[0]) {
289d25d7f95SJed Brown #if defined(PETSC_SIZEOF_LONG_LONG)
290d25d7f95SJed Brown         long long lid;
291d25d7f95SJed Brown         sscanf(str,"%lld",&lid);
2921a1499c8SBarry 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);
293c272547aSJed Brown #else
294d25d7f95SJed Brown         long  lid;
295d25d7f95SJed Brown         sscanf(str,"%ld",&lid);
2961a1499c8SBarry 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);
297c272547aSJed Brown #endif
298a297a907SKarl Rupp 
299d25d7f95SJed Brown         next->set = PETSC_TRUE;
300d25d7f95SJed Brown         *((PetscInt*)next->data) = (PetscInt)lid;
301aee2cecaSBarry Smith       }
302aee2cecaSBarry Smith       break;
303aee2cecaSBarry Smith     case OPTION_REAL:
304e55864a3SBarry 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);
3053fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3063fc1eb6aSBarry Smith       if (str[0]) {
307ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
308a4404d99SBarry Smith         sscanf(str,"%e",&ir);
309570b7f6dSBarry Smith #elif defined(PETSC_USE_REAL___FP16)
310570b7f6dSBarry Smith         float irtemp;
311570b7f6dSBarry Smith         sscanf(str,"%e",&irtemp);
312570b7f6dSBarry Smith         ir = irtemp;
313ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL_DOUBLE)
314aee2cecaSBarry Smith         sscanf(str,"%le",&ir);
315ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL___FLOAT128)
316d9822059SBarry Smith         ir = strtoflt128(str,0);
317d9822059SBarry Smith #else
318513dbe71SLisandro Dalcin         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Unknown scalar type");
319a4404d99SBarry Smith #endif
320aee2cecaSBarry Smith         next->set                 = PETSC_TRUE;
321aee2cecaSBarry Smith         *((PetscReal*)next->data) = ir;
322aee2cecaSBarry Smith       }
323aee2cecaSBarry Smith       break;
3247781c08eSBarry Smith     case OPTION_BOOL:
32583355fc5SBarry 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);
3267781c08eSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3277781c08eSBarry Smith       if (str[0]) {
3287781c08eSBarry Smith         ierr = PetscOptionsStringToBool(str,&bid);CHKERRQ(ierr);
3297781c08eSBarry Smith         next->set = PETSC_TRUE;
3307781c08eSBarry Smith         *((PetscBool*)next->data) = bid;
3317781c08eSBarry Smith       }
3327781c08eSBarry Smith       break;
333aee2cecaSBarry Smith     case OPTION_STRING:
334e55864a3SBarry 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);
3353fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3363fc1eb6aSBarry Smith       if (str[0]) {
337aee2cecaSBarry Smith         next->set = PETSC_TRUE;
33864facd6cSBarry Smith         /* must use system malloc since SAWs may free this */
3395b02f95dSBarry Smith         ierr = PetscStrdup(str,(char**)&next->data);CHKERRQ(ierr);
3406356e834SBarry Smith       }
3416356e834SBarry Smith       break;
342a264d7a6SBarry Smith     case OPTION_FLIST:
343e55864a3SBarry Smith       ierr = PetscFunctionListPrintTypes(PETSC_COMM_WORLD,stdout,PetscOptionsObject->prefix,next->option,next->text,next->man,next->flist,(char*)next->data);CHKERRQ(ierr);
3443cc1e11dSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3453cc1e11dSBarry Smith       if (str[0]) {
346e55864a3SBarry Smith         PetscOptionsObject->changedmethod = PETSC_TRUE;
3473cc1e11dSBarry Smith         next->set = PETSC_TRUE;
34864facd6cSBarry Smith         /* must use system malloc since SAWs may free this */
3495b02f95dSBarry Smith         ierr = PetscStrdup(str,(char**)&next->data);CHKERRQ(ierr);
3503cc1e11dSBarry Smith       }
3513cc1e11dSBarry Smith       break;
352b432afdaSMatthew Knepley     default:
353b432afdaSMatthew Knepley       break;
3546356e834SBarry Smith     }
3556356e834SBarry Smith     next = next->next;
3566356e834SBarry Smith   }
3576356e834SBarry Smith   PetscFunctionReturn(0);
3586356e834SBarry Smith }
3596356e834SBarry Smith 
360e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
361e04113cfSBarry Smith #include <petscviewersaws.h>
362d5649816SBarry Smith 
363d5649816SBarry Smith static int count = 0;
364d5649816SBarry Smith 
365e04113cfSBarry Smith PetscErrorCode PetscOptionsSAWsDestroy(void)
366d5649816SBarry Smith {
3672657e9d9SBarry Smith   PetscFunctionBegin;
368d5649816SBarry Smith   PetscFunctionReturn(0);
369d5649816SBarry Smith }
370d5649816SBarry Smith 
3719c1e0ce8SBarry Smith static const char *OptionsHeader = "<head>\n"
37223a1ff79SBarry Smith                                    "<script type=\"text/javascript\" src=\"http://www.mcs.anl.gov/research/projects/saws/js/jquery-1.9.1.js\"></script>\n"
37323a1ff79SBarry Smith                                    "<script type=\"text/javascript\" src=\"http://www.mcs.anl.gov/research/projects/saws/js/SAWs.js\"></script>\n"
374d1fc0251SBarry Smith                                    "<script type=\"text/javascript\" src=\"js/PETSc.js\"></script>\n"
37564bbc9efSBarry Smith                                    "<script>\n"
37664bbc9efSBarry Smith                                       "jQuery(document).ready(function() {\n"
37764bbc9efSBarry Smith                                          "PETSc.getAndDisplayDirectory(null,\"#variablesInfo\")\n"
37864bbc9efSBarry Smith                                       "})\n"
37964bbc9efSBarry Smith                                   "</script>\n"
38064bbc9efSBarry Smith                                   "</head>\n";
3811423471aSBarry Smith 
3821423471aSBarry Smith /*  Determines the size and style of the scroll region where PETSc options selectable from users are displayed */
3831423471aSBarry Smith static const char *OptionsBodyBottom = "<div id=\"variablesInfo\" style=\"background-color:lightblue;height:auto;max-height:500px;overflow:scroll;\"></div>\n<br>\n</body>";
38464bbc9efSBarry Smith 
385b3506946SBarry Smith /*
3867781c08eSBarry Smith     PetscOptionsSAWsInput - Presents all the PETSc Options processed by the program so the user may change them at runtime using the SAWs
387b3506946SBarry Smith 
388b3506946SBarry Smith     Bugs:
389b3506946SBarry Smith +    All processes must traverse through the exact same set of option queries do to the call to PetscScanString()
390b3506946SBarry Smith .    Internal strings have arbitrary length and string copies are not checked that they fit into string space
391b3506946SBarry Smith -    Only works for PetscInt == int, PetscReal == double etc
392b3506946SBarry Smith 
393b3506946SBarry Smith 
394b3506946SBarry Smith */
3954416b707SBarry Smith PetscErrorCode PetscOptionsSAWsInput(PetscOptionItems *PetscOptionsObject)
396b3506946SBarry Smith {
397b3506946SBarry Smith   PetscErrorCode  ierr;
3984416b707SBarry Smith   PetscOptionItem next     = PetscOptionsObject->next;
399d5649816SBarry Smith   static int      mancount = 0;
400b3506946SBarry Smith   char            options[16];
4017aab2a10SBarry Smith   PetscBool       changedmethod = PETSC_FALSE;
402a23eb982SSurtai Han   PetscBool       stopasking    = PETSC_FALSE;
40388a9752cSBarry Smith   char            manname[16],textname[16];
4042657e9d9SBarry Smith   char            dir[1024];
405b3506946SBarry Smith 
4062a409bb0SBarry Smith   PetscFunctionBegin;
407b3506946SBarry Smith   /* the next line is a bug, this will only work if all processors are here, the comm passed in is ignored!!! */
408b3506946SBarry Smith   sprintf(options,"Options_%d",count++);
409a297a907SKarl Rupp 
4107325c4a4SBarry Smith   PetscOptionsObject->pprefix = PetscOptionsObject->prefix; /* SAWs will change this, so cannot pass prefix directly */
4111bc75a8dSBarry Smith 
412d50831c4SBarry Smith   ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s","_title");CHKERRQ(ierr);
41383355fc5SBarry Smith   PetscStackCallSAWs(SAWs_Register,(dir,&PetscOptionsObject->title,1,SAWs_READ,SAWs_STRING));
4147781c08eSBarry Smith   ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s","prefix");CHKERRQ(ierr);
41583355fc5SBarry Smith   PetscStackCallSAWs(SAWs_Register,(dir,&PetscOptionsObject->pprefix,1,SAWs_READ,SAWs_STRING));
4162657e9d9SBarry Smith   PetscStackCallSAWs(SAWs_Register,("/PETSc/Options/ChangedMethod",&changedmethod,1,SAWs_WRITE,SAWs_BOOLEAN));
417a23eb982SSurtai Han   PetscStackCallSAWs(SAWs_Register,("/PETSc/Options/StopAsking",&stopasking,1,SAWs_WRITE,SAWs_BOOLEAN));
418b3506946SBarry Smith 
419b3506946SBarry Smith   while (next) {
420d50831c4SBarry Smith     sprintf(manname,"_man_%d",mancount);
4212657e9d9SBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",manname);CHKERRQ(ierr);
4227781c08eSBarry Smith     PetscStackCallSAWs(SAWs_Register,(dir,&next->man,1,SAWs_READ,SAWs_STRING));
423d50831c4SBarry Smith     sprintf(textname,"_text_%d",mancount++);
4247781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",textname);CHKERRQ(ierr);
4257781c08eSBarry Smith     PetscStackCallSAWs(SAWs_Register,(dir,&next->text,1,SAWs_READ,SAWs_STRING));
4269f32e415SBarry Smith 
427b3506946SBarry Smith     switch (next->type) {
428b3506946SBarry Smith     case OPTION_HEAD:
429b3506946SBarry Smith       break;
430b3506946SBarry Smith     case OPTION_INT_ARRAY:
4317781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4322657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_INT));
433b3506946SBarry Smith       break;
434b3506946SBarry Smith     case OPTION_REAL_ARRAY:
4357781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4362657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_DOUBLE));
437b3506946SBarry Smith       break;
438b3506946SBarry Smith     case OPTION_INT:
4397781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4402657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_INT));
441b3506946SBarry Smith       break;
442b3506946SBarry Smith     case OPTION_REAL:
4437781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4442657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_DOUBLE));
445b3506946SBarry Smith       break;
4467781c08eSBarry Smith     case OPTION_BOOL:
4477781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4482657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_BOOLEAN));
4491ae3d29cSBarry Smith       break;
4507781c08eSBarry Smith     case OPTION_BOOL_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_BOOLEAN));
45371f08665SBarry Smith       break;
454b3506946SBarry Smith     case OPTION_STRING:
4557781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4567781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
4571ae3d29cSBarry Smith       break;
4581ae3d29cSBarry Smith     case OPTION_STRING_ARRAY:
4597781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4602657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_STRING));
461b3506946SBarry Smith       break;
462a264d7a6SBarry Smith     case OPTION_FLIST:
463a264d7a6SBarry Smith       {
464a264d7a6SBarry Smith       PetscInt ntext;
4657781c08eSBarry Smith       ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4667781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
467a264d7a6SBarry Smith       ierr = PetscFunctionListGet(next->flist,(const char***)&next->edata,&ntext);CHKERRQ(ierr);
468a264d7a6SBarry Smith       PetscStackCallSAWs(SAWs_Set_Legal_Variable_Values,(dir,ntext,next->edata));
469a264d7a6SBarry Smith       }
470a264d7a6SBarry Smith       break;
4711ae3d29cSBarry Smith     case OPTION_ELIST:
472a264d7a6SBarry Smith       {
473a264d7a6SBarry Smith       PetscInt ntext = next->nlist;
4747781c08eSBarry Smith       ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4757781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
476ead66b60SBarry Smith       ierr = PetscMalloc1((ntext+1),(char***)&next->edata);CHKERRQ(ierr);
4771ae3d29cSBarry Smith       ierr = PetscMemcpy(next->edata,next->list,ntext*sizeof(char*));CHKERRQ(ierr);
478a264d7a6SBarry Smith       PetscStackCallSAWs(SAWs_Set_Legal_Variable_Values,(dir,ntext,next->edata));
479a264d7a6SBarry Smith       }
480a264d7a6SBarry Smith       break;
481b3506946SBarry Smith     default:
482b3506946SBarry Smith       break;
483b3506946SBarry Smith     }
484b3506946SBarry Smith     next = next->next;
485b3506946SBarry Smith   }
486b3506946SBarry Smith 
487b3506946SBarry Smith   /* wait until accessor has unlocked the memory */
48864bbc9efSBarry Smith   PetscStackCallSAWs(SAWs_Push_Header,("index.html",OptionsHeader));
48964bbc9efSBarry Smith   PetscStackCallSAWs(SAWs_Push_Body,("index.html",2,OptionsBodyBottom));
4907aab2a10SBarry Smith   ierr = PetscSAWsBlock();CHKERRQ(ierr);
49164bbc9efSBarry Smith   PetscStackCallSAWs(SAWs_Pop_Header,("index.html"));
49264bbc9efSBarry Smith   PetscStackCallSAWs(SAWs_Pop_Body,("index.html",2));
493b3506946SBarry Smith 
49488a9752cSBarry Smith   /* determine if any values have been set in GUI */
49583355fc5SBarry Smith   next = PetscOptionsObject->next;
49688a9752cSBarry Smith   while (next) {
49788a9752cSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
498f7b25cf6SBarry Smith     PetscStackCallSAWs(SAWs_Selected,(dir,(int*)&next->set));
49988a9752cSBarry Smith     next = next->next;
50088a9752cSBarry Smith   }
50188a9752cSBarry Smith 
502b3506946SBarry Smith   /* reset counter to -2; this updates the screen with the new options for the selected method */
503f7b25cf6SBarry Smith   if (changedmethod) PetscOptionsObject->count = -2;
504b3506946SBarry Smith 
505a23eb982SSurtai Han   if (stopasking) {
506a23eb982SSurtai Han     PetscOptionsPublish      = PETSC_FALSE;
50712655325SBarry Smith     PetscOptionsObject->count = 0;//do not ask for same thing again
508a23eb982SSurtai Han   }
509a23eb982SSurtai Han 
5109a492a5cSBarry Smith   PetscStackCallSAWs(SAWs_Delete,("/PETSc/Options"));
511b3506946SBarry Smith   PetscFunctionReturn(0);
512b3506946SBarry Smith }
513b3506946SBarry Smith #endif
514b3506946SBarry Smith 
5154416b707SBarry Smith PetscErrorCode PetscOptionsEnd_Private(PetscOptionItems *PetscOptionsObject)
51653acd3b1SBarry Smith {
51753acd3b1SBarry Smith   PetscErrorCode  ierr;
5184416b707SBarry Smith   PetscOptionItem last;
5196356e834SBarry Smith   char            option[256],value[1024],tmp[32];
520330cf3c9SBarry Smith   size_t          j;
52153acd3b1SBarry Smith 
52253acd3b1SBarry Smith   PetscFunctionBegin;
52383355fc5SBarry Smith   if (PetscOptionsObject->next) {
52483355fc5SBarry Smith     if (!PetscOptionsObject->count) {
525a264d7a6SBarry Smith #if defined(PETSC_HAVE_SAWS)
526f7b25cf6SBarry Smith       ierr = PetscOptionsSAWsInput(PetscOptionsObject);CHKERRQ(ierr);
527b3506946SBarry Smith #else
528e55864a3SBarry Smith       ierr = PetscOptionsGetFromTextInput(PetscOptionsObject);CHKERRQ(ierr);
529b3506946SBarry Smith #endif
530aee2cecaSBarry Smith     }
531aee2cecaSBarry Smith   }
5326356e834SBarry Smith 
533e55864a3SBarry Smith   ierr = PetscFree(PetscOptionsObject->title);CHKERRQ(ierr);
5346356e834SBarry Smith 
535e26ddf31SBarry Smith   /* reset counter to -2; this updates the screen with the new options for the selected method */
536e55864a3SBarry Smith   if (PetscOptionsObject->changedmethod) PetscOptionsObject->count = -2;
5377a72a596SBarry Smith   /* reset alreadyprinted flag */
538e55864a3SBarry Smith   PetscOptionsObject->alreadyprinted = PETSC_FALSE;
539e55864a3SBarry Smith   if (PetscOptionsObject->object) PetscOptionsObject->object->optionsprinted = PETSC_TRUE;
540e55864a3SBarry Smith   PetscOptionsObject->object = NULL;
54153acd3b1SBarry Smith 
542e55864a3SBarry Smith   while (PetscOptionsObject->next) {
543e55864a3SBarry Smith     if (PetscOptionsObject->next->set) {
544e55864a3SBarry Smith       if (PetscOptionsObject->prefix) {
54553acd3b1SBarry Smith         ierr = PetscStrcpy(option,"-");CHKERRQ(ierr);
546e55864a3SBarry Smith         ierr = PetscStrcat(option,PetscOptionsObject->prefix);CHKERRQ(ierr);
547e55864a3SBarry Smith         ierr = PetscStrcat(option,PetscOptionsObject->next->option+1);CHKERRQ(ierr);
5486356e834SBarry Smith       } else {
549e55864a3SBarry Smith         ierr = PetscStrcpy(option,PetscOptionsObject->next->option);CHKERRQ(ierr);
5506356e834SBarry Smith       }
5516356e834SBarry Smith 
552e55864a3SBarry Smith       switch (PetscOptionsObject->next->type) {
5536356e834SBarry Smith       case OPTION_HEAD:
5546356e834SBarry Smith         break;
5556356e834SBarry Smith       case OPTION_INT_ARRAY:
556e55864a3SBarry Smith         sprintf(value,"%d",(int)((PetscInt*)PetscOptionsObject->next->data)[0]);
557e55864a3SBarry Smith         for (j=1; j<PetscOptionsObject->next->arraylength; j++) {
558e55864a3SBarry Smith           sprintf(tmp,"%d",(int)((PetscInt*)PetscOptionsObject->next->data)[j]);
5596356e834SBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5606356e834SBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5616356e834SBarry Smith         }
5626356e834SBarry Smith         break;
5636356e834SBarry Smith       case OPTION_INT:
564e55864a3SBarry Smith         sprintf(value,"%d",(int) *(PetscInt*)PetscOptionsObject->next->data);
5656356e834SBarry Smith         break;
5666356e834SBarry Smith       case OPTION_REAL:
567e55864a3SBarry Smith         sprintf(value,"%g",(double) *(PetscReal*)PetscOptionsObject->next->data);
5686356e834SBarry Smith         break;
5696356e834SBarry Smith       case OPTION_REAL_ARRAY:
570e55864a3SBarry Smith         sprintf(value,"%g",(double)((PetscReal*)PetscOptionsObject->next->data)[0]);
571e55864a3SBarry Smith         for (j=1; j<PetscOptionsObject->next->arraylength; j++) {
572e55864a3SBarry Smith           sprintf(tmp,"%g",(double)((PetscReal*)PetscOptionsObject->next->data)[j]);
5736356e834SBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5746356e834SBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5756356e834SBarry Smith         }
5766356e834SBarry Smith         break;
577050cccc3SHong Zhang       case OPTION_SCALAR_ARRAY:
57895f3a755SJose E. Roman         sprintf(value,"%g+%gi",(double)PetscRealPart(((PetscScalar*)PetscOptionsObject->next->data)[0]),(double)PetscImaginaryPart(((PetscScalar*)PetscOptionsObject->next->data)[0]));
579050cccc3SHong Zhang         for (j=1; j<PetscOptionsObject->next->arraylength; j++) {
58095f3a755SJose E. Roman           sprintf(tmp,"%g+%gi",(double)PetscRealPart(((PetscScalar*)PetscOptionsObject->next->data)[j]),(double)PetscImaginaryPart(((PetscScalar*)PetscOptionsObject->next->data)[j]));
581050cccc3SHong Zhang           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
582050cccc3SHong Zhang           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
583050cccc3SHong Zhang         }
584050cccc3SHong Zhang         break;
5857781c08eSBarry Smith       case OPTION_BOOL:
586e55864a3SBarry Smith         sprintf(value,"%d",*(int*)PetscOptionsObject->next->data);
5876356e834SBarry Smith         break;
5887781c08eSBarry Smith       case OPTION_BOOL_ARRAY:
589e55864a3SBarry Smith         sprintf(value,"%d",(int)((PetscBool*)PetscOptionsObject->next->data)[0]);
590e55864a3SBarry Smith         for (j=1; j<PetscOptionsObject->next->arraylength; j++) {
591e55864a3SBarry Smith           sprintf(tmp,"%d",(int)((PetscBool*)PetscOptionsObject->next->data)[j]);
5921ae3d29cSBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5931ae3d29cSBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5941ae3d29cSBarry Smith         }
5951ae3d29cSBarry Smith         break;
596a264d7a6SBarry Smith       case OPTION_FLIST:
5976991f827SBarry Smith         ierr = PetscStrcpy(value,(char*)PetscOptionsObject->next->data);CHKERRQ(ierr);
5986991f827SBarry Smith         break;
5991ae3d29cSBarry Smith       case OPTION_ELIST:
600e55864a3SBarry Smith         ierr = PetscStrcpy(value,(char*)PetscOptionsObject->next->data);CHKERRQ(ierr);
6016356e834SBarry Smith         break;
6021ae3d29cSBarry Smith       case OPTION_STRING:
603e55864a3SBarry Smith         ierr = PetscStrcpy(value,(char*)PetscOptionsObject->next->data);CHKERRQ(ierr);
604d51da6bfSBarry Smith         break;
6051ae3d29cSBarry Smith       case OPTION_STRING_ARRAY:
606e55864a3SBarry Smith         sprintf(value,"%s",((char**)PetscOptionsObject->next->data)[0]);
607e55864a3SBarry Smith         for (j=1; j<PetscOptionsObject->next->arraylength; j++) {
608e55864a3SBarry Smith           sprintf(tmp,"%s",((char**)PetscOptionsObject->next->data)[j]);
6091ae3d29cSBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
6101ae3d29cSBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
6111ae3d29cSBarry Smith         }
6126356e834SBarry Smith         break;
6136356e834SBarry Smith       }
614c5929fdfSBarry Smith       ierr = PetscOptionsSetValue(PetscOptionsObject->options,option,value);CHKERRQ(ierr);
6156356e834SBarry Smith     }
6166991f827SBarry Smith     if (PetscOptionsObject->next->type == OPTION_ELIST) {
6176991f827SBarry Smith       ierr = PetscStrNArrayDestroy(PetscOptionsObject->next->nlist,(char ***)&PetscOptionsObject->next->list);CHKERRQ(ierr);
6186991f827SBarry Smith     }
619e55864a3SBarry Smith     ierr   = PetscFree(PetscOptionsObject->next->text);CHKERRQ(ierr);
620e55864a3SBarry Smith     ierr   = PetscFree(PetscOptionsObject->next->option);CHKERRQ(ierr);
621e55864a3SBarry Smith     ierr   = PetscFree(PetscOptionsObject->next->man);CHKERRQ(ierr);
622e55864a3SBarry Smith     ierr   = PetscFree(PetscOptionsObject->next->edata);CHKERRQ(ierr);
623c979a496SBarry Smith 
62483355fc5SBarry Smith     if ((PetscOptionsObject->next->type == OPTION_STRING) || (PetscOptionsObject->next->type == OPTION_FLIST) || (PetscOptionsObject->next->type == OPTION_ELIST)){
62583355fc5SBarry Smith       free(PetscOptionsObject->next->data);
626c979a496SBarry Smith     } else {
62783355fc5SBarry Smith       ierr   = PetscFree(PetscOptionsObject->next->data);CHKERRQ(ierr);
628c979a496SBarry Smith     }
6297781c08eSBarry Smith 
63083355fc5SBarry Smith     last                    = PetscOptionsObject->next;
63183355fc5SBarry Smith     PetscOptionsObject->next = PetscOptionsObject->next->next;
6326356e834SBarry Smith     ierr                    = PetscFree(last);CHKERRQ(ierr);
6336356e834SBarry Smith   }
634f59f755dSBarry Smith   ierr = PetscFree(PetscOptionsObject->prefix);CHKERRQ(ierr);
635e55864a3SBarry Smith   PetscOptionsObject->next = 0;
63653acd3b1SBarry Smith   PetscFunctionReturn(0);
63753acd3b1SBarry Smith }
63853acd3b1SBarry Smith 
63953acd3b1SBarry Smith /*@C
64053acd3b1SBarry Smith    PetscOptionsEnum - Gets the enum value for a particular option in the database.
64153acd3b1SBarry Smith 
6423f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
64353acd3b1SBarry Smith 
64453acd3b1SBarry Smith    Input Parameters:
64553acd3b1SBarry Smith +  opt - option name
64653acd3b1SBarry Smith .  text - short string that describes the option
64753acd3b1SBarry Smith .  man - manual page with additional information on option
64853acd3b1SBarry Smith .  list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null
6490fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
6500fdccdaeSBarry Smith $                 PetscOptionsEnum(..., obj->value,&object->value,...) or
6510fdccdaeSBarry Smith $                 value = defaultvalue
6520fdccdaeSBarry Smith $                 PetscOptionsEnum(..., value,&value,&flg);
6530fdccdaeSBarry Smith $                 if (flg) {
65453acd3b1SBarry Smith 
65553acd3b1SBarry Smith    Output Parameter:
65653acd3b1SBarry Smith +  value - the  value to return
657b32e0204SMatthew G Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
65853acd3b1SBarry Smith 
65953acd3b1SBarry Smith    Level: beginner
66053acd3b1SBarry Smith 
66153acd3b1SBarry Smith    Concepts: options database
66253acd3b1SBarry Smith 
66353acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
66453acd3b1SBarry Smith 
66553acd3b1SBarry Smith           list is usually something like PCASMTypes or some other predefined list of enum names
66653acd3b1SBarry Smith 
6672efd9cb1SBarry Smith           If the user does not supply the option at all value is NOT changed. Thus
6682efd9cb1SBarry Smith           you should ALWAYS initialize value if you access it without first checking if the set flag is true.
6692efd9cb1SBarry Smith 
670*989712b9SBarry Smith           The default/currentvalue passed into this routine does not get transferred to the output value variable automatically.
671*989712b9SBarry Smith 
67253acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
673acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
674acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
67553acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
67653acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
677acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
678a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
67953acd3b1SBarry Smith @*/
6804416b707SBarry 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)
68153acd3b1SBarry Smith {
68253acd3b1SBarry Smith   PetscErrorCode ierr;
68353acd3b1SBarry Smith   PetscInt       ntext = 0;
684aa5bb8c0SSatish Balay   PetscInt       tval;
685ace3abfcSBarry Smith   PetscBool      tflg;
68653acd3b1SBarry Smith 
68753acd3b1SBarry Smith   PetscFunctionBegin;
68853acd3b1SBarry Smith   while (list[ntext++]) {
689e32f2f54SBarry Smith     if (ntext > 50) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument appears to be wrong or have more than 50 entries");
69053acd3b1SBarry Smith   }
691e32f2f54SBarry Smith   if (ntext < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument must have at least two entries: typename and type prefix");
69253acd3b1SBarry Smith   ntext -= 3;
693e55864a3SBarry Smith   ierr   = PetscOptionsEList_Private(PetscOptionsObject,opt,text,man,list,ntext,list[currentvalue],&tval,&tflg);CHKERRQ(ierr);
694aa5bb8c0SSatish Balay   /* with PETSC_USE_64BIT_INDICES sizeof(PetscInt) != sizeof(PetscEnum) */
695aa5bb8c0SSatish Balay   if (tflg) *value = (PetscEnum)tval;
696aa5bb8c0SSatish Balay   if (set)  *set   = tflg;
69753acd3b1SBarry Smith   PetscFunctionReturn(0);
69853acd3b1SBarry Smith }
69953acd3b1SBarry Smith 
700d3e47460SLisandro Dalcin /*@C
701d3e47460SLisandro Dalcin    PetscOptionsEnumArray - Gets an array of enum values for a particular
702d3e47460SLisandro Dalcin    option in the database.
703d3e47460SLisandro Dalcin 
704d3e47460SLisandro Dalcin    Logically Collective on the communicator passed in PetscOptionsBegin()
705d3e47460SLisandro Dalcin 
706d3e47460SLisandro Dalcin    Input Parameters:
707d3e47460SLisandro Dalcin +  opt - the option one is seeking
708d3e47460SLisandro Dalcin .  text - short string describing option
709d3e47460SLisandro Dalcin .  man - manual page for option
710d3e47460SLisandro Dalcin -  n - maximum number of values
711d3e47460SLisandro Dalcin 
712d3e47460SLisandro Dalcin    Output Parameter:
713d3e47460SLisandro Dalcin +  value - location to copy values
714d3e47460SLisandro Dalcin .  n - actual number of values found
715d3e47460SLisandro Dalcin -  set - PETSC_TRUE if found, else PETSC_FALSE
716d3e47460SLisandro Dalcin 
717d3e47460SLisandro Dalcin    Level: beginner
718d3e47460SLisandro Dalcin 
719d3e47460SLisandro Dalcin    Notes:
720d3e47460SLisandro Dalcin    The array must be passed as a comma separated list.
721d3e47460SLisandro Dalcin 
722d3e47460SLisandro Dalcin    There must be no intervening spaces between the values.
723d3e47460SLisandro Dalcin 
724d3e47460SLisandro Dalcin    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
725d3e47460SLisandro Dalcin 
726d3e47460SLisandro Dalcin    Concepts: options database^array of enums
727d3e47460SLisandro Dalcin 
728d3e47460SLisandro Dalcin .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
729d3e47460SLisandro Dalcin           PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
730d3e47460SLisandro Dalcin           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
731d3e47460SLisandro Dalcin           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
732d3e47460SLisandro Dalcin           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
733d3e47460SLisandro Dalcin           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsRealArray()
734d3e47460SLisandro Dalcin @*/
7354416b707SBarry 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)
736d3e47460SLisandro Dalcin {
737d3e47460SLisandro Dalcin   PetscInt        i,nlist = 0;
7384416b707SBarry Smith   PetscOptionItem amsopt;
739d3e47460SLisandro Dalcin   PetscErrorCode  ierr;
740d3e47460SLisandro Dalcin 
741d3e47460SLisandro Dalcin   PetscFunctionBegin;
742d3e47460SLisandro 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");
743d3e47460SLisandro Dalcin   if (nlist < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument must have at least two entries: typename and type prefix");
744d3e47460SLisandro Dalcin   nlist -= 3; /* drop enum name, prefix, and null termination */
745d3e47460SLisandro Dalcin   if (0 && !PetscOptionsObject->count) { /* XXX Requires additional support */
746d3e47460SLisandro Dalcin     PetscEnum *vals;
7474416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_INT_ARRAY/*XXX OPTION_ENUM_ARRAY*/,&amsopt);CHKERRQ(ierr);
748d3e47460SLisandro Dalcin     ierr = PetscStrNArrayallocpy(nlist,list,(char***)&amsopt->list);CHKERRQ(ierr);
749d3e47460SLisandro Dalcin     amsopt->nlist = nlist;
750d3e47460SLisandro Dalcin     ierr = PetscMalloc1(*n,(PetscEnum**)&amsopt->data);CHKERRQ(ierr);
751d3e47460SLisandro Dalcin     amsopt->arraylength = *n;
752d3e47460SLisandro Dalcin     vals = (PetscEnum*)amsopt->data;
753d3e47460SLisandro Dalcin     for (i=0; i<*n; i++) vals[i] = value[i];
754d3e47460SLisandro Dalcin   }
755c5929fdfSBarry Smith   ierr = PetscOptionsGetEnumArray(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,list,value,n,set);CHKERRQ(ierr);
756d3e47460SLisandro Dalcin   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
757d3e47460SLisandro Dalcin     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%s",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,list[value[0]]);CHKERRQ(ierr);
758d3e47460SLisandro Dalcin     for (i=1; i<*n; i++) {ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%s",list[value[i]]);CHKERRQ(ierr);}
759d3e47460SLisandro Dalcin     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (choose from)",text);CHKERRQ(ierr);
760d3e47460SLisandro Dalcin     for (i=0; i<nlist; i++) {ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," %s",list[i]);CHKERRQ(ierr);}
761d3e47460SLisandro Dalcin     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," (%s)\n",ManSection(man));CHKERRQ(ierr);
762d3e47460SLisandro Dalcin   }
763d3e47460SLisandro Dalcin   PetscFunctionReturn(0);
764d3e47460SLisandro Dalcin }
765d3e47460SLisandro Dalcin 
76653acd3b1SBarry Smith /* -------------------------------------------------------------------------------------------------------------*/
76753acd3b1SBarry Smith /*@C
76853acd3b1SBarry Smith    PetscOptionsInt - Gets the integer value for a particular option in the database.
76953acd3b1SBarry Smith 
7703f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
77153acd3b1SBarry Smith 
77253acd3b1SBarry Smith    Input Parameters:
77353acd3b1SBarry Smith +  opt - option name
77453acd3b1SBarry Smith .  text - short string that describes the option
77553acd3b1SBarry Smith .  man - manual page with additional information on option
7760fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
7770fdccdaeSBarry Smith $                 PetscOptionsInt(..., obj->value,&object->value,...) or
7780fdccdaeSBarry Smith $                 value = defaultvalue
7790fdccdaeSBarry Smith $                 PetscOptionsInt(..., value,&value,&flg);
7800fdccdaeSBarry Smith $                 if (flg) {
78153acd3b1SBarry Smith 
78253acd3b1SBarry Smith    Output Parameter:
78353acd3b1SBarry Smith +  value - the integer value to return
78453acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
78553acd3b1SBarry Smith 
7862efd9cb1SBarry Smith    Notes: If the user does not supply the option at all value is NOT changed. Thus
7872efd9cb1SBarry Smith           you should ALWAYS initialize value if you access it without first checking if the set flag is true.
7882efd9cb1SBarry Smith 
789*989712b9SBarry Smith           The default/currentvalue passed into this routine does not get transferred to the output value variable automatically.
790*989712b9SBarry Smith 
79153acd3b1SBarry Smith    Level: beginner
79253acd3b1SBarry Smith 
79353acd3b1SBarry Smith    Concepts: options database^has int
79453acd3b1SBarry Smith 
79553acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
79653acd3b1SBarry Smith 
79753acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
798acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
799acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
80053acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
80153acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
802acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
803a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
80453acd3b1SBarry Smith @*/
8054416b707SBarry Smith PetscErrorCode  PetscOptionsInt_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscInt currentvalue,PetscInt *value,PetscBool  *set)
80653acd3b1SBarry Smith {
80753acd3b1SBarry Smith   PetscErrorCode  ierr;
8084416b707SBarry Smith   PetscOptionItem amsopt;
80912655325SBarry Smith   PetscBool       wasset;
81053acd3b1SBarry Smith 
81153acd3b1SBarry Smith   PetscFunctionBegin;
812e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
8134416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_INT,&amsopt);CHKERRQ(ierr);
8146356e834SBarry Smith     ierr = PetscMalloc(sizeof(PetscInt),&amsopt->data);CHKERRQ(ierr);
81512655325SBarry Smith     *(PetscInt*)amsopt->data = currentvalue;
8163e211508SBarry Smith 
817c5929fdfSBarry Smith     ierr = PetscOptionsGetInt(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,&currentvalue,&wasset);CHKERRQ(ierr);
8183e211508SBarry Smith     if (wasset) {
81912655325SBarry Smith       *(PetscInt*)amsopt->data = currentvalue;
8203e211508SBarry Smith     }
821af6d86caSBarry Smith   }
822c5929fdfSBarry Smith   ierr = PetscOptionsGetInt(NULL,PetscOptionsObject->prefix,opt,value,set);CHKERRQ(ierr);
823e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
8241a1499c8SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%d>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,currentvalue,text,ManSection(man));CHKERRQ(ierr);
82553acd3b1SBarry Smith   }
82653acd3b1SBarry Smith   PetscFunctionReturn(0);
82753acd3b1SBarry Smith }
82853acd3b1SBarry Smith 
82953acd3b1SBarry Smith /*@C
83053acd3b1SBarry Smith    PetscOptionsString - Gets the string value for a particular option in the database.
83153acd3b1SBarry Smith 
8323f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
83353acd3b1SBarry Smith 
83453acd3b1SBarry Smith    Input Parameters:
83553acd3b1SBarry Smith +  opt - option name
83653acd3b1SBarry Smith .  text - short string that describes the option
83753acd3b1SBarry Smith .  man - manual page with additional information on option
8380fdccdaeSBarry Smith .  currentvalue - the current value; caller is responsible for setting this value correctly. This is not used to set value
839bcbf2dc5SJed Brown -  len - length of the result string including null terminator
84053acd3b1SBarry Smith 
84153acd3b1SBarry Smith    Output Parameter:
84253acd3b1SBarry Smith +  value - the value to return
84353acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
84453acd3b1SBarry Smith 
84553acd3b1SBarry Smith    Level: beginner
84653acd3b1SBarry Smith 
84753acd3b1SBarry Smith    Concepts: options database^has int
84853acd3b1SBarry Smith 
84953acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
85053acd3b1SBarry Smith 
8517fccdfe4SBarry 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).
8527fccdfe4SBarry Smith 
8532efd9cb1SBarry Smith           If the user does not supply the option at all value is NOT changed. Thus
8542efd9cb1SBarry Smith           you should ALWAYS initialize value if you access it without first checking if the set flag is true.
8552efd9cb1SBarry Smith 
856*989712b9SBarry Smith           The default/currentvalue passed into this routine does not get transferred to the output value variable automatically.
857*989712b9SBarry Smith 
8582efd9cb1SBarry Smith 
859c5929fdfSBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(NULL,),
860acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
861acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsReal(), PetscOptionsBool(),
86253acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
86353acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
864acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
865a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
86653acd3b1SBarry Smith @*/
8674416b707SBarry 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)
86853acd3b1SBarry Smith {
86953acd3b1SBarry Smith   PetscErrorCode  ierr;
8704416b707SBarry Smith   PetscOptionItem amsopt;
87153acd3b1SBarry Smith 
87253acd3b1SBarry Smith   PetscFunctionBegin;
8731a1499c8SBarry Smith   if (!PetscOptionsObject->count) {
8744416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
87564facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
8760fdccdaeSBarry Smith     ierr = PetscStrdup(currentvalue ? currentvalue : "",(char**)&amsopt->data);CHKERRQ(ierr);
877af6d86caSBarry Smith   }
878c5929fdfSBarry Smith   ierr = PetscOptionsGetString(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,len,set);CHKERRQ(ierr);
879e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
8801a1499c8SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%s>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,currentvalue,text,ManSection(man));CHKERRQ(ierr);
88153acd3b1SBarry Smith   }
88253acd3b1SBarry Smith   PetscFunctionReturn(0);
88353acd3b1SBarry Smith }
88453acd3b1SBarry Smith 
88553acd3b1SBarry Smith /*@C
88653acd3b1SBarry Smith    PetscOptionsReal - Gets the PetscReal value for a particular option in the database.
88753acd3b1SBarry Smith 
8883f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
88953acd3b1SBarry Smith 
89053acd3b1SBarry Smith    Input Parameters:
89153acd3b1SBarry Smith +  opt - option name
89253acd3b1SBarry Smith .  text - short string that describes the option
89353acd3b1SBarry Smith .  man - manual page with additional information on option
8940fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
8950fdccdaeSBarry Smith $                 PetscOptionsReal(..., obj->value,&object->value,...) or
8960fdccdaeSBarry Smith $                 value = defaultvalue
8970fdccdaeSBarry Smith $                 PetscOptionsReal(..., value,&value,&flg);
8980fdccdaeSBarry Smith $                 if (flg) {
89953acd3b1SBarry Smith 
90053acd3b1SBarry Smith    Output Parameter:
90153acd3b1SBarry Smith +  value - the value to return
90253acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
90353acd3b1SBarry Smith 
9042efd9cb1SBarry Smith    Notes:  If the user does not supply the option at all value is NOT changed. Thus
9052efd9cb1SBarry Smith           you should ALWAYS initialize value if you access it without first checking if the set flag is true.
9062efd9cb1SBarry Smith 
907*989712b9SBarry Smith           The default/currentvalue passed into this routine does not get transferred to the output value variable automatically.
908*989712b9SBarry Smith 
90953acd3b1SBarry Smith    Level: beginner
91053acd3b1SBarry Smith 
91153acd3b1SBarry Smith    Concepts: options database^has int
91253acd3b1SBarry Smith 
91353acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
91453acd3b1SBarry Smith 
915c5929fdfSBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(NULL,),
916acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
917acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
91853acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
91953acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
920acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
921a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
92253acd3b1SBarry Smith @*/
9234416b707SBarry Smith PetscErrorCode  PetscOptionsReal_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscReal currentvalue,PetscReal *value,PetscBool  *set)
92453acd3b1SBarry Smith {
92553acd3b1SBarry Smith   PetscErrorCode  ierr;
9264416b707SBarry Smith   PetscOptionItem amsopt;
92753acd3b1SBarry Smith 
92853acd3b1SBarry Smith   PetscFunctionBegin;
929e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
9304416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_REAL,&amsopt);CHKERRQ(ierr);
931538aa990SBarry Smith     ierr = PetscMalloc(sizeof(PetscReal),&amsopt->data);CHKERRQ(ierr);
932a297a907SKarl Rupp 
9330fdccdaeSBarry Smith     *(PetscReal*)amsopt->data = currentvalue;
934538aa990SBarry Smith   }
935c5929fdfSBarry Smith   ierr = PetscOptionsGetReal(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,set);CHKERRQ(ierr);
9361a1499c8SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
9371a1499c8SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%g>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,(double)currentvalue,text,ManSection(man));CHKERRQ(ierr);
93853acd3b1SBarry Smith   }
93953acd3b1SBarry Smith   PetscFunctionReturn(0);
94053acd3b1SBarry Smith }
94153acd3b1SBarry Smith 
94253acd3b1SBarry Smith /*@C
94353acd3b1SBarry Smith    PetscOptionsScalar - Gets the scalar value for a particular option in the database.
94453acd3b1SBarry Smith 
9453f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
94653acd3b1SBarry Smith 
94753acd3b1SBarry Smith    Input Parameters:
94853acd3b1SBarry Smith +  opt - option name
94953acd3b1SBarry Smith .  text - short string that describes the option
95053acd3b1SBarry Smith .  man - manual page with additional information on option
9510fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
9520fdccdaeSBarry Smith $                 PetscOptionsScalar(..., obj->value,&object->value,...) or
9530fdccdaeSBarry Smith $                 value = defaultvalue
9540fdccdaeSBarry Smith $                 PetscOptionsScalar(..., value,&value,&flg);
9550fdccdaeSBarry Smith $                 if (flg) {
9560fdccdaeSBarry Smith 
95753acd3b1SBarry Smith 
95853acd3b1SBarry Smith    Output Parameter:
95953acd3b1SBarry Smith +  value - the value to return
96053acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
96153acd3b1SBarry Smith 
9622efd9cb1SBarry Smith    Notes: If the user does not supply the option at all value is NOT changed. Thus
9632efd9cb1SBarry Smith           you should ALWAYS initialize value if you access it without first checking if the set flag is true.
9642efd9cb1SBarry Smith 
965*989712b9SBarry Smith           The default/currentvalue passed into this routine does not get transferred to the output value variable automatically.
966*989712b9SBarry Smith 
96753acd3b1SBarry Smith    Level: beginner
96853acd3b1SBarry Smith 
96953acd3b1SBarry Smith    Concepts: options database^has int
97053acd3b1SBarry Smith 
97153acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
97253acd3b1SBarry Smith 
973c5929fdfSBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(NULL,),
974acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
975acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
97653acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
97753acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
978acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
979a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
98053acd3b1SBarry Smith @*/
9814416b707SBarry Smith PetscErrorCode  PetscOptionsScalar_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscScalar currentvalue,PetscScalar *value,PetscBool  *set)
98253acd3b1SBarry Smith {
98353acd3b1SBarry Smith   PetscErrorCode ierr;
98453acd3b1SBarry Smith 
98553acd3b1SBarry Smith   PetscFunctionBegin;
98653acd3b1SBarry Smith #if !defined(PETSC_USE_COMPLEX)
9870fdccdaeSBarry Smith   ierr = PetscOptionsReal(opt,text,man,currentvalue,value,set);CHKERRQ(ierr);
98853acd3b1SBarry Smith #else
989c5929fdfSBarry Smith   ierr = PetscOptionsGetScalar(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,set);CHKERRQ(ierr);
99053acd3b1SBarry Smith #endif
99153acd3b1SBarry Smith   PetscFunctionReturn(0);
99253acd3b1SBarry Smith }
99353acd3b1SBarry Smith 
99453acd3b1SBarry Smith /*@C
99590d69ab7SBarry Smith    PetscOptionsName - Determines if a particular option has been set in the database. This returns true whether the option is a number, string or boolean, even
99690d69ab7SBarry Smith                       its value is set to false.
99753acd3b1SBarry Smith 
9983f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
99953acd3b1SBarry Smith 
100053acd3b1SBarry Smith    Input Parameters:
100153acd3b1SBarry Smith +  opt - option name
100253acd3b1SBarry Smith .  text - short string that describes the option
100353acd3b1SBarry Smith -  man - manual page with additional information on option
100453acd3b1SBarry Smith 
100553acd3b1SBarry Smith    Output Parameter:
100653acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
100753acd3b1SBarry Smith 
100853acd3b1SBarry Smith    Level: beginner
100953acd3b1SBarry Smith 
101053acd3b1SBarry Smith    Concepts: options database^has int
101153acd3b1SBarry Smith 
101253acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
101353acd3b1SBarry Smith 
1014c5929fdfSBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(NULL,),
1015acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
1016acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
101753acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
101853acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1019acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1020a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
102153acd3b1SBarry Smith @*/
10224416b707SBarry Smith PetscErrorCode  PetscOptionsName_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool  *flg)
102353acd3b1SBarry Smith {
102453acd3b1SBarry Smith   PetscErrorCode  ierr;
10254416b707SBarry Smith   PetscOptionItem amsopt;
102653acd3b1SBarry Smith 
102753acd3b1SBarry Smith   PetscFunctionBegin;
1028e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
10294416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1030ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1031a297a907SKarl Rupp 
1032ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
10331ae3d29cSBarry Smith   }
1034c5929fdfSBarry Smith   ierr = PetscOptionsHasName(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,flg);CHKERRQ(ierr);
1035e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1036e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
103753acd3b1SBarry Smith   }
103853acd3b1SBarry Smith   PetscFunctionReturn(0);
103953acd3b1SBarry Smith }
104053acd3b1SBarry Smith 
104153acd3b1SBarry Smith /*@C
1042a264d7a6SBarry Smith      PetscOptionsFList - Puts a list of option values that a single one may be selected from
104353acd3b1SBarry Smith 
10443f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
104553acd3b1SBarry Smith 
104653acd3b1SBarry Smith    Input Parameters:
104753acd3b1SBarry Smith +  opt - option name
104853acd3b1SBarry Smith .  text - short string that describes the option
104953acd3b1SBarry Smith .  man - manual page with additional information on option
105053acd3b1SBarry Smith .  list - the possible choices
10510fdccdaeSBarry Smith .  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with
10520fdccdaeSBarry Smith $                 PetscOptionsFlist(..., obj->value,value,len,&flg);
10530fdccdaeSBarry Smith $                 if (flg) {
10543cc1e11dSBarry Smith -  len - the length of the character array value
105553acd3b1SBarry Smith 
105653acd3b1SBarry Smith    Output Parameter:
105753acd3b1SBarry Smith +  value - the value to return
105853acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
105953acd3b1SBarry Smith 
106053acd3b1SBarry Smith    Level: intermediate
106153acd3b1SBarry Smith 
106253acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
106353acd3b1SBarry Smith 
10642efd9cb1SBarry Smith           If the user does not supply the option at all value is NOT changed. Thus
10652efd9cb1SBarry Smith           you should ALWAYS initialize value if you access it without first checking if the set flag is true.
10662efd9cb1SBarry Smith 
1067*989712b9SBarry Smith           The default/currentvalue passed into this routine does not get transferred to the output value variable automatically.
1068*989712b9SBarry Smith 
106953acd3b1SBarry Smith    See PetscOptionsEList() for when the choices are given in a string array
107053acd3b1SBarry Smith 
107153acd3b1SBarry Smith    To get a listing of all currently specified options,
107288c29154SBarry Smith     see PetscOptionsView() or PetscOptionsGetAll()
107353acd3b1SBarry Smith 
1074eabe10d7SBarry Smith    Developer Note: This cannot check for invalid selection because of things like MATAIJ that are not included in the list
1075eabe10d7SBarry Smith 
107653acd3b1SBarry Smith    Concepts: options database^list
107753acd3b1SBarry Smith 
1078c5929fdfSBarry Smith .seealso: PetscOptionsGetInt(NULL,), PetscOptionsGetReal(),
1079acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
108053acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
108153acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1082acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1083a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsEnum()
108453acd3b1SBarry Smith @*/
10854416b707SBarry 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)
108653acd3b1SBarry Smith {
108753acd3b1SBarry Smith   PetscErrorCode  ierr;
10884416b707SBarry Smith   PetscOptionItem amsopt;
108953acd3b1SBarry Smith 
109053acd3b1SBarry Smith   PetscFunctionBegin;
10911a1499c8SBarry Smith   if (!PetscOptionsObject->count) {
10924416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,ltext,man,OPTION_FLIST,&amsopt);CHKERRQ(ierr);
109364facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
10940fdccdaeSBarry Smith     ierr = PetscStrdup(currentvalue ? currentvalue : "",(char**)&amsopt->data);CHKERRQ(ierr);
10953cc1e11dSBarry Smith     amsopt->flist = list;
10963cc1e11dSBarry Smith   }
1097c5929fdfSBarry Smith   ierr = PetscOptionsGetString(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,len,set);CHKERRQ(ierr);
10981a1499c8SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
10991a1499c8SBarry Smith     ierr = PetscFunctionListPrintTypes(PetscOptionsObject->comm,stdout,PetscOptionsObject->prefix,opt,ltext,man,list,currentvalue);CHKERRQ(ierr);CHKERRQ(ierr);
110053acd3b1SBarry Smith   }
110153acd3b1SBarry Smith   PetscFunctionReturn(0);
110253acd3b1SBarry Smith }
110353acd3b1SBarry Smith 
110453acd3b1SBarry Smith /*@C
110553acd3b1SBarry Smith      PetscOptionsEList - Puts a list of option values that a single one may be selected from
110653acd3b1SBarry Smith 
11073f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
110853acd3b1SBarry Smith 
110953acd3b1SBarry Smith    Input Parameters:
111053acd3b1SBarry Smith +  opt - option name
111153acd3b1SBarry Smith .  ltext - short string that describes the option
111253acd3b1SBarry Smith .  man - manual page with additional information on option
1113a264d7a6SBarry Smith .  list - the possible choices (one of these must be selected, anything else is invalid)
111453acd3b1SBarry Smith .  ntext - number of choices
11150fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with
11160fdccdaeSBarry Smith $                 PetscOptionsElist(..., obj->value,&value,&flg);
11170fdccdaeSBarry Smith $                 if (flg) {
11180fdccdaeSBarry Smith 
111953acd3b1SBarry Smith 
112053acd3b1SBarry Smith    Output Parameter:
112153acd3b1SBarry Smith +  value - the index of the value to return
112253acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
112353acd3b1SBarry Smith 
112453acd3b1SBarry Smith    Level: intermediate
112553acd3b1SBarry Smith 
112653acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
112753acd3b1SBarry Smith 
11282efd9cb1SBarry Smith          If the user does not supply the option at all value is NOT changed. Thus
11292efd9cb1SBarry Smith           you should ALWAYS initialize value if you access it without first checking if the set flag is true.
11302efd9cb1SBarry Smith 
1131a264d7a6SBarry Smith    See PetscOptionsFList() for when the choices are given in a PetscFunctionList()
113253acd3b1SBarry Smith 
113353acd3b1SBarry Smith    Concepts: options database^list
113453acd3b1SBarry Smith 
1135c5929fdfSBarry Smith .seealso: PetscOptionsGetInt(NULL,), PetscOptionsGetReal(),
1136acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
113753acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
113853acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1139acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1140a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEnum()
114153acd3b1SBarry Smith @*/
11424416b707SBarry 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)
114353acd3b1SBarry Smith {
114453acd3b1SBarry Smith   PetscErrorCode  ierr;
114553acd3b1SBarry Smith   PetscInt        i;
11464416b707SBarry Smith   PetscOptionItem amsopt;
114753acd3b1SBarry Smith 
114853acd3b1SBarry Smith   PetscFunctionBegin;
11491a1499c8SBarry Smith   if (!PetscOptionsObject->count) {
11504416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,ltext,man,OPTION_ELIST,&amsopt);CHKERRQ(ierr);
115164facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
11520fdccdaeSBarry Smith     ierr = PetscStrdup(currentvalue ? currentvalue : "",(char**)&amsopt->data);CHKERRQ(ierr);
11536991f827SBarry Smith     ierr = PetscStrNArrayallocpy(ntext,list,(char***)&amsopt->list);CHKERRQ(ierr);
11541ae3d29cSBarry Smith     amsopt->nlist = ntext;
11551ae3d29cSBarry Smith   }
1156c5929fdfSBarry Smith   ierr = PetscOptionsGetEList(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,list,ntext,value,set);CHKERRQ(ierr);
11571a1499c8SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1158e3f729a5SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%s> %s (choose one of)",PetscOptionsObject->prefix?PetscOptionsObject->prefix:"",opt+1,currentvalue,ltext);CHKERRQ(ierr);
115953acd3b1SBarry Smith     for (i=0; i<ntext; i++) {
1160e55864a3SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," %s",list[i]);CHKERRQ(ierr);
116153acd3b1SBarry Smith     }
1162e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," (%s)\n",ManSection(man));CHKERRQ(ierr);
116353acd3b1SBarry Smith   }
116453acd3b1SBarry Smith   PetscFunctionReturn(0);
116553acd3b1SBarry Smith }
116653acd3b1SBarry Smith 
116753acd3b1SBarry Smith /*@C
1168acfcf0e5SJed Brown      PetscOptionsBoolGroupBegin - First in a series of logical queries on the options database for
1169d5649816SBarry Smith        which at most a single value can be true.
117053acd3b1SBarry Smith 
11713f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
117253acd3b1SBarry Smith 
117353acd3b1SBarry Smith    Input Parameters:
117453acd3b1SBarry Smith +  opt - option name
117553acd3b1SBarry Smith .  text - short string that describes the option
117653acd3b1SBarry Smith -  man - manual page with additional information on option
117753acd3b1SBarry Smith 
117853acd3b1SBarry Smith    Output Parameter:
117953acd3b1SBarry Smith .  flg - whether that option was set or not
118053acd3b1SBarry Smith 
118153acd3b1SBarry Smith    Level: intermediate
118253acd3b1SBarry Smith 
118353acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
118453acd3b1SBarry Smith 
1185acfcf0e5SJed Brown    Must be followed by 0 or more PetscOptionsBoolGroup()s and PetscOptionsBoolGroupEnd()
118653acd3b1SBarry Smith 
118753acd3b1SBarry Smith     Concepts: options database^logical group
118853acd3b1SBarry Smith 
1189c5929fdfSBarry Smith .seealso: PetscOptionsGetInt(NULL,), PetscOptionsGetReal(),
1190acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
119153acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
119253acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1193acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1194a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
119553acd3b1SBarry Smith @*/
11964416b707SBarry Smith PetscErrorCode  PetscOptionsBoolGroupBegin_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool  *flg)
119753acd3b1SBarry Smith {
119853acd3b1SBarry Smith   PetscErrorCode  ierr;
11994416b707SBarry Smith   PetscOptionItem amsopt;
120053acd3b1SBarry Smith 
120153acd3b1SBarry Smith   PetscFunctionBegin;
1202e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
12034416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1204ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1205a297a907SKarl Rupp 
1206ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
12071ae3d29cSBarry Smith   }
120868b16fdaSBarry Smith   *flg = PETSC_FALSE;
1209c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,flg,NULL);CHKERRQ(ierr);
1210e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1211e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  Pick at most one of -------------\n");CHKERRQ(ierr);
1212e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"    -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
121353acd3b1SBarry Smith   }
121453acd3b1SBarry Smith   PetscFunctionReturn(0);
121553acd3b1SBarry Smith }
121653acd3b1SBarry Smith 
121753acd3b1SBarry Smith /*@C
1218acfcf0e5SJed Brown      PetscOptionsBoolGroup - One in a series of logical queries on the options database for
1219d5649816SBarry Smith        which at most a single value can be true.
122053acd3b1SBarry Smith 
12213f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
122253acd3b1SBarry Smith 
122353acd3b1SBarry Smith    Input Parameters:
122453acd3b1SBarry Smith +  opt - option name
122553acd3b1SBarry Smith .  text - short string that describes the option
122653acd3b1SBarry Smith -  man - manual page with additional information on option
122753acd3b1SBarry Smith 
122853acd3b1SBarry Smith    Output Parameter:
122953acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
123053acd3b1SBarry Smith 
123153acd3b1SBarry Smith    Level: intermediate
123253acd3b1SBarry Smith 
123353acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
123453acd3b1SBarry Smith 
1235acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin() and preceded a PetscOptionsBoolGroupEnd()
123653acd3b1SBarry Smith 
123753acd3b1SBarry Smith     Concepts: options database^logical group
123853acd3b1SBarry Smith 
1239c5929fdfSBarry Smith .seealso: PetscOptionsGetInt(NULL,), PetscOptionsGetReal(),
1240acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
124153acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
124253acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1243acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1244a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
124553acd3b1SBarry Smith @*/
12464416b707SBarry Smith PetscErrorCode  PetscOptionsBoolGroup_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool  *flg)
124753acd3b1SBarry Smith {
124853acd3b1SBarry Smith   PetscErrorCode  ierr;
12494416b707SBarry Smith   PetscOptionItem amsopt;
125053acd3b1SBarry Smith 
125153acd3b1SBarry Smith   PetscFunctionBegin;
1252e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
12534416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1254ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1255a297a907SKarl Rupp 
1256ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
12571ae3d29cSBarry Smith   }
125817326d04SJed Brown   *flg = PETSC_FALSE;
1259c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,flg,NULL);CHKERRQ(ierr);
1260e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1261e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"    -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
126253acd3b1SBarry Smith   }
126353acd3b1SBarry Smith   PetscFunctionReturn(0);
126453acd3b1SBarry Smith }
126553acd3b1SBarry Smith 
126653acd3b1SBarry Smith /*@C
1267acfcf0e5SJed Brown      PetscOptionsBoolGroupEnd - Last in a series of logical queries on the options database for
1268d5649816SBarry Smith        which at most a single value can be true.
126953acd3b1SBarry Smith 
12703f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
127153acd3b1SBarry Smith 
127253acd3b1SBarry Smith    Input Parameters:
127353acd3b1SBarry Smith +  opt - option name
127453acd3b1SBarry Smith .  text - short string that describes the option
127553acd3b1SBarry Smith -  man - manual page with additional information on option
127653acd3b1SBarry Smith 
127753acd3b1SBarry Smith    Output Parameter:
127853acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
127953acd3b1SBarry Smith 
128053acd3b1SBarry Smith    Level: intermediate
128153acd3b1SBarry Smith 
128253acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
128353acd3b1SBarry Smith 
1284acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin()
128553acd3b1SBarry Smith 
128653acd3b1SBarry Smith     Concepts: options database^logical group
128753acd3b1SBarry Smith 
1288c5929fdfSBarry Smith .seealso: PetscOptionsGetInt(NULL,), PetscOptionsGetReal(),
1289acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
129053acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
129153acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1292acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1293a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
129453acd3b1SBarry Smith @*/
12954416b707SBarry Smith PetscErrorCode  PetscOptionsBoolGroupEnd_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool  *flg)
129653acd3b1SBarry Smith {
129753acd3b1SBarry Smith   PetscErrorCode  ierr;
12984416b707SBarry Smith   PetscOptionItem amsopt;
129953acd3b1SBarry Smith 
130053acd3b1SBarry Smith   PetscFunctionBegin;
1301e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
13024416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1303ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1304a297a907SKarl Rupp 
1305ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
13061ae3d29cSBarry Smith   }
130717326d04SJed Brown   *flg = PETSC_FALSE;
1308c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,flg,NULL);CHKERRQ(ierr);
1309e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1310e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"    -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
131153acd3b1SBarry Smith   }
131253acd3b1SBarry Smith   PetscFunctionReturn(0);
131353acd3b1SBarry Smith }
131453acd3b1SBarry Smith 
131553acd3b1SBarry Smith /*@C
1316acfcf0e5SJed Brown    PetscOptionsBool - Determines if a particular option is in the database with a true or false
131753acd3b1SBarry Smith 
13183f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
131953acd3b1SBarry Smith 
132053acd3b1SBarry Smith    Input Parameters:
132153acd3b1SBarry Smith +  opt - option name
132253acd3b1SBarry Smith .  text - short string that describes the option
1323868c398cSBarry Smith .  man - manual page with additional information on option
132494ae4db5SBarry Smith -  currentvalue - the current value
132553acd3b1SBarry Smith 
132653acd3b1SBarry Smith    Output Parameter:
132753acd3b1SBarry Smith .  flg - PETSC_TRUE or PETSC_FALSE
132853acd3b1SBarry Smith .  set - PETSC_TRUE if found, else PETSC_FALSE
132953acd3b1SBarry Smith 
13302efd9cb1SBarry Smith    Notes:
13312efd9cb1SBarry Smith        TRUE, true, YES, yes, nostring, and 1 all translate to PETSC_TRUE
13322efd9cb1SBarry Smith        FALSE, false, NO, no, and 0 all translate to PETSC_FALSE
13332efd9cb1SBarry Smith 
13342efd9cb1SBarry Smith       If the option is given, but no value is provided, then flg and set are both given the value PETSC_TRUE. That is -requested_bool
13352efd9cb1SBarry Smith      is equivalent to -requested_bool true
13362efd9cb1SBarry Smith 
13372efd9cb1SBarry Smith        If the user does not supply the option at all flg is NOT changed. Thus
13382efd9cb1SBarry Smith      you should ALWAYS initialize the flg if you access it without first checking if the set flag is true.
13392efd9cb1SBarry 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 /*@C
137753acd3b1SBarry Smith    PetscOptionsRealArray - Gets an array of double values for a particular
137853acd3b1SBarry Smith    option in the database. The values must be separated with commas with
137953acd3b1SBarry Smith    no intervening spaces.
138053acd3b1SBarry Smith 
13813f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
138253acd3b1SBarry Smith 
138353acd3b1SBarry Smith    Input Parameters:
138453acd3b1SBarry Smith +  opt - the option one is seeking
138553acd3b1SBarry Smith .  text - short string describing option
138653acd3b1SBarry Smith .  man - manual page for option
138753acd3b1SBarry Smith -  nmax - maximum number of values
138853acd3b1SBarry Smith 
138953acd3b1SBarry Smith    Output Parameter:
139053acd3b1SBarry Smith +  value - location to copy values
139153acd3b1SBarry Smith .  nmax - actual number of values found
139253acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
139353acd3b1SBarry Smith 
139453acd3b1SBarry Smith    Level: beginner
139553acd3b1SBarry Smith 
139653acd3b1SBarry Smith    Notes:
139753acd3b1SBarry Smith    The user should pass in an array of doubles
139853acd3b1SBarry Smith 
139953acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
140053acd3b1SBarry Smith 
140153acd3b1SBarry Smith    Concepts: options database^array of strings
140253acd3b1SBarry Smith 
1403c5929fdfSBarry Smith .seealso: PetscOptionsGetInt(NULL,), PetscOptionsGetReal(),
1404acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
140553acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
140653acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1407acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1408a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
140953acd3b1SBarry Smith @*/
14104416b707SBarry Smith PetscErrorCode PetscOptionsRealArray_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscReal value[],PetscInt *n,PetscBool  *set)
141153acd3b1SBarry Smith {
141253acd3b1SBarry Smith   PetscErrorCode  ierr;
141353acd3b1SBarry Smith   PetscInt        i;
14144416b707SBarry Smith   PetscOptionItem amsopt;
141553acd3b1SBarry Smith 
141653acd3b1SBarry Smith   PetscFunctionBegin;
1417e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
1418e26ddf31SBarry Smith     PetscReal *vals;
1419e26ddf31SBarry Smith 
14204416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_REAL_ARRAY,&amsopt);CHKERRQ(ierr);
1421e55864a3SBarry Smith     ierr = PetscMalloc((*n)*sizeof(PetscReal),&amsopt->data);CHKERRQ(ierr);
1422e26ddf31SBarry Smith     vals = (PetscReal*)amsopt->data;
1423e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1424e26ddf31SBarry Smith     amsopt->arraylength = *n;
1425e26ddf31SBarry Smith   }
1426c5929fdfSBarry Smith   ierr = PetscOptionsGetRealArray(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,n,set);CHKERRQ(ierr);
1427e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1428a519f713SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%g",PetscOptionsObject->prefix?PetscOptionsObject->prefix:"",opt+1,(double)value[0]);CHKERRQ(ierr);
142953acd3b1SBarry Smith     for (i=1; i<*n; i++) {
1430e55864a3SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%g",(double)value[i]);CHKERRQ(ierr);
143153acd3b1SBarry Smith     }
1432e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
143353acd3b1SBarry Smith   }
143453acd3b1SBarry Smith   PetscFunctionReturn(0);
143553acd3b1SBarry Smith }
143653acd3b1SBarry Smith 
1437050cccc3SHong Zhang /*@C
1438050cccc3SHong Zhang    PetscOptionsScalarArray - Gets an array of Scalar values for a particular
1439050cccc3SHong Zhang    option in the database. The values must be separated with commas with
1440050cccc3SHong Zhang    no intervening spaces.
1441050cccc3SHong Zhang 
1442050cccc3SHong Zhang    Logically Collective on the communicator passed in PetscOptionsBegin()
1443050cccc3SHong Zhang 
1444050cccc3SHong Zhang    Input Parameters:
1445050cccc3SHong Zhang +  opt - the option one is seeking
1446050cccc3SHong Zhang .  text - short string describing option
1447050cccc3SHong Zhang .  man - manual page for option
1448050cccc3SHong Zhang -  nmax - maximum number of values
1449050cccc3SHong Zhang 
1450050cccc3SHong Zhang    Output Parameter:
1451050cccc3SHong Zhang +  value - location to copy values
1452050cccc3SHong Zhang .  nmax - actual number of values found
1453050cccc3SHong Zhang -  set - PETSC_TRUE if found, else PETSC_FALSE
1454050cccc3SHong Zhang 
1455050cccc3SHong Zhang    Level: beginner
1456050cccc3SHong Zhang 
1457050cccc3SHong Zhang    Notes:
1458050cccc3SHong Zhang    The user should pass in an array of doubles
1459050cccc3SHong Zhang 
1460050cccc3SHong Zhang    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
1461050cccc3SHong Zhang 
1462050cccc3SHong Zhang    Concepts: options database^array of strings
1463050cccc3SHong Zhang 
1464c5929fdfSBarry Smith .seealso: PetscOptionsGetInt(NULL,), PetscOptionsGetReal(),
1465050cccc3SHong Zhang           PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
1466050cccc3SHong Zhang           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1467050cccc3SHong Zhang           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1468050cccc3SHong Zhang           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1469050cccc3SHong Zhang           PetscOptionsFList(), PetscOptionsEList()
1470050cccc3SHong Zhang @*/
14714416b707SBarry Smith PetscErrorCode PetscOptionsScalarArray_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscScalar value[],PetscInt *n,PetscBool  *set)
1472050cccc3SHong Zhang {
1473050cccc3SHong Zhang   PetscErrorCode  ierr;
1474050cccc3SHong Zhang   PetscInt        i;
14754416b707SBarry Smith   PetscOptionItem amsopt;
1476050cccc3SHong Zhang 
1477050cccc3SHong Zhang   PetscFunctionBegin;
1478050cccc3SHong Zhang   if (!PetscOptionsObject->count) {
1479050cccc3SHong Zhang     PetscScalar *vals;
1480050cccc3SHong Zhang 
14814416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_SCALAR_ARRAY,&amsopt);CHKERRQ(ierr);
1482050cccc3SHong Zhang     ierr = PetscMalloc((*n)*sizeof(PetscScalar),&amsopt->data);CHKERRQ(ierr);
1483050cccc3SHong Zhang     vals = (PetscScalar*)amsopt->data;
1484050cccc3SHong Zhang     for (i=0; i<*n; i++) vals[i] = value[i];
1485050cccc3SHong Zhang     amsopt->arraylength = *n;
1486050cccc3SHong Zhang   }
1487c5929fdfSBarry Smith   ierr = PetscOptionsGetScalarArray(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,n,set);CHKERRQ(ierr);
1488050cccc3SHong Zhang   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
148995f3a755SJose 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);
1490050cccc3SHong Zhang     for (i=1; i<*n; i++) {
149195f3a755SJose E. Roman       ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%g+%gi",(double)PetscRealPart(value[i]),(double)PetscImaginaryPart(value[i]));CHKERRQ(ierr);
1492050cccc3SHong Zhang     }
1493050cccc3SHong Zhang     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
1494050cccc3SHong Zhang   }
1495050cccc3SHong Zhang   PetscFunctionReturn(0);
1496050cccc3SHong Zhang }
149753acd3b1SBarry Smith 
149853acd3b1SBarry Smith /*@C
149953acd3b1SBarry Smith    PetscOptionsIntArray - Gets an array of integers for a particular
1500b32a342fSShri Abhyankar    option in the database.
150153acd3b1SBarry Smith 
15023f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
150353acd3b1SBarry Smith 
150453acd3b1SBarry Smith    Input Parameters:
150553acd3b1SBarry Smith +  opt - the option one is seeking
150653acd3b1SBarry Smith .  text - short string describing option
150753acd3b1SBarry Smith .  man - manual page for option
1508f8a50e2bSBarry Smith -  n - maximum number of values
150953acd3b1SBarry Smith 
151053acd3b1SBarry Smith    Output Parameter:
151153acd3b1SBarry Smith +  value - location to copy values
1512f8a50e2bSBarry Smith .  n - actual number of values found
151353acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
151453acd3b1SBarry Smith 
151553acd3b1SBarry Smith    Level: beginner
151653acd3b1SBarry Smith 
151753acd3b1SBarry Smith    Notes:
1518b32a342fSShri Abhyankar    The array can be passed as
1519bebe2cf6SSatish Balay    a comma separated list:                                 0,1,2,3,4,5,6,7
15200fd488f5SShri Abhyankar    a range (start-end+1):                                  0-8
15210fd488f5SShri Abhyankar    a range with given increment (start-end+1:inc):         0-7:2
1522bebe2cf6SSatish Balay    a combination of values and ranges separated by commas: 0,1-8,8-15:2
1523b32a342fSShri Abhyankar 
1524b32a342fSShri Abhyankar    There must be no intervening spaces between the values.
152553acd3b1SBarry Smith 
152653acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
152753acd3b1SBarry Smith 
1528b32a342fSShri Abhyankar    Concepts: options database^array of ints
152953acd3b1SBarry Smith 
1530c5929fdfSBarry Smith .seealso: PetscOptionsGetInt(NULL,), PetscOptionsGetReal(),
1531acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
153253acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
153353acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1534acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1535a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsRealArray()
153653acd3b1SBarry Smith @*/
15374416b707SBarry Smith PetscErrorCode  PetscOptionsIntArray_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscInt value[],PetscInt *n,PetscBool  *set)
153853acd3b1SBarry Smith {
153953acd3b1SBarry Smith   PetscErrorCode ierr;
154053acd3b1SBarry Smith   PetscInt        i;
15414416b707SBarry Smith   PetscOptionItem amsopt;
154253acd3b1SBarry Smith 
154353acd3b1SBarry Smith   PetscFunctionBegin;
1544e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
1545e26ddf31SBarry Smith     PetscInt *vals;
1546e26ddf31SBarry Smith 
15474416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_INT_ARRAY,&amsopt);CHKERRQ(ierr);
1548854ce69bSBarry Smith     ierr = PetscMalloc1(*n,(PetscInt**)&amsopt->data);CHKERRQ(ierr);
1549e26ddf31SBarry Smith     vals = (PetscInt*)amsopt->data;
1550e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1551e26ddf31SBarry Smith     amsopt->arraylength = *n;
1552e26ddf31SBarry Smith   }
1553c5929fdfSBarry Smith   ierr = PetscOptionsGetIntArray(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,n,set);CHKERRQ(ierr);
1554e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1555e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%d",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,value[0]);CHKERRQ(ierr);
155653acd3b1SBarry Smith     for (i=1; i<*n; i++) {
1557e55864a3SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%d",value[i]);CHKERRQ(ierr);
155853acd3b1SBarry Smith     }
1559e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
156053acd3b1SBarry Smith   }
156153acd3b1SBarry Smith   PetscFunctionReturn(0);
156253acd3b1SBarry Smith }
156353acd3b1SBarry Smith 
156453acd3b1SBarry Smith /*@C
156553acd3b1SBarry Smith    PetscOptionsStringArray - Gets an array of string values for a particular
156653acd3b1SBarry Smith    option in the database. The values must be separated with commas with
156753acd3b1SBarry Smith    no intervening spaces.
156853acd3b1SBarry Smith 
15693f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
157053acd3b1SBarry Smith 
157153acd3b1SBarry Smith    Input Parameters:
157253acd3b1SBarry Smith +  opt - the option one is seeking
157353acd3b1SBarry Smith .  text - short string describing option
157453acd3b1SBarry Smith .  man - manual page for option
157553acd3b1SBarry Smith -  nmax - maximum number of strings
157653acd3b1SBarry Smith 
157753acd3b1SBarry Smith    Output Parameter:
157853acd3b1SBarry Smith +  value - location to copy strings
157953acd3b1SBarry Smith .  nmax - actual number of strings found
158053acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
158153acd3b1SBarry Smith 
158253acd3b1SBarry Smith    Level: beginner
158353acd3b1SBarry Smith 
158453acd3b1SBarry Smith    Notes:
158553acd3b1SBarry Smith    The user should pass in an array of pointers to char, to hold all the
158653acd3b1SBarry Smith    strings returned by this function.
158753acd3b1SBarry Smith 
158853acd3b1SBarry Smith    The user is responsible for deallocating the strings that are
158953acd3b1SBarry Smith    returned. The Fortran interface for this routine is not supported.
159053acd3b1SBarry Smith 
159153acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
159253acd3b1SBarry Smith 
159353acd3b1SBarry Smith    Concepts: options database^array of strings
159453acd3b1SBarry Smith 
1595c5929fdfSBarry Smith .seealso: PetscOptionsGetInt(NULL,), PetscOptionsGetReal(),
1596acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
159753acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
159853acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1599acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1600a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
160153acd3b1SBarry Smith @*/
16024416b707SBarry Smith PetscErrorCode  PetscOptionsStringArray_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],char *value[],PetscInt *nmax,PetscBool  *set)
160353acd3b1SBarry Smith {
160453acd3b1SBarry Smith   PetscErrorCode  ierr;
16054416b707SBarry Smith   PetscOptionItem amsopt;
160653acd3b1SBarry Smith 
160753acd3b1SBarry Smith   PetscFunctionBegin;
1608e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
16094416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_STRING_ARRAY,&amsopt);CHKERRQ(ierr);
1610854ce69bSBarry Smith     ierr = PetscMalloc1(*nmax,(char**)&amsopt->data);CHKERRQ(ierr);
1611a297a907SKarl Rupp 
16121ae3d29cSBarry Smith     amsopt->arraylength = *nmax;
16131ae3d29cSBarry Smith   }
1614c5929fdfSBarry Smith   ierr = PetscOptionsGetStringArray(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,nmax,set);CHKERRQ(ierr);
1615e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1616e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <string1,string2,...>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
161753acd3b1SBarry Smith   }
161853acd3b1SBarry Smith   PetscFunctionReturn(0);
161953acd3b1SBarry Smith }
162053acd3b1SBarry Smith 
1621e2446a98SMatthew Knepley /*@C
1622acfcf0e5SJed Brown    PetscOptionsBoolArray - Gets an array of logical values (true or false) for a particular
1623e2446a98SMatthew Knepley    option in the database. The values must be separated with commas with
1624e2446a98SMatthew Knepley    no intervening spaces.
1625e2446a98SMatthew Knepley 
16263f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
1627e2446a98SMatthew Knepley 
1628e2446a98SMatthew Knepley    Input Parameters:
1629e2446a98SMatthew Knepley +  opt - the option one is seeking
1630e2446a98SMatthew Knepley .  text - short string describing option
1631e2446a98SMatthew Knepley .  man - manual page for option
1632e2446a98SMatthew Knepley -  nmax - maximum number of values
1633e2446a98SMatthew Knepley 
1634e2446a98SMatthew Knepley    Output Parameter:
1635e2446a98SMatthew Knepley +  value - location to copy values
1636e2446a98SMatthew Knepley .  nmax - actual number of values found
1637e2446a98SMatthew Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
1638e2446a98SMatthew Knepley 
1639e2446a98SMatthew Knepley    Level: beginner
1640e2446a98SMatthew Knepley 
1641e2446a98SMatthew Knepley    Notes:
1642e2446a98SMatthew Knepley    The user should pass in an array of doubles
1643e2446a98SMatthew Knepley 
1644e2446a98SMatthew Knepley    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
1645e2446a98SMatthew Knepley 
1646e2446a98SMatthew Knepley    Concepts: options database^array of strings
1647e2446a98SMatthew Knepley 
1648c5929fdfSBarry Smith .seealso: PetscOptionsGetInt(NULL,), PetscOptionsGetReal(),
1649acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
1650e2446a98SMatthew Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1651e2446a98SMatthew Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1652acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1653a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
1654e2446a98SMatthew Knepley @*/
16554416b707SBarry Smith PetscErrorCode  PetscOptionsBoolArray_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool value[],PetscInt *n,PetscBool *set)
1656e2446a98SMatthew Knepley {
1657e2446a98SMatthew Knepley   PetscErrorCode   ierr;
1658e2446a98SMatthew Knepley   PetscInt         i;
16594416b707SBarry Smith   PetscOptionItem  amsopt;
1660e2446a98SMatthew Knepley 
1661e2446a98SMatthew Knepley   PetscFunctionBegin;
1662e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
1663ace3abfcSBarry Smith     PetscBool *vals;
16641ae3d29cSBarry Smith 
16654416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL_ARRAY,&amsopt);CHKERRQ(ierr);
16661a1499c8SBarry Smith     ierr = PetscMalloc1(*n,(PetscBool**)&amsopt->data);CHKERRQ(ierr);
1667ace3abfcSBarry Smith     vals = (PetscBool*)amsopt->data;
16681ae3d29cSBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
16691ae3d29cSBarry Smith     amsopt->arraylength = *n;
16701ae3d29cSBarry Smith   }
1671c5929fdfSBarry Smith   ierr = PetscOptionsGetBoolArray(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,n,set);CHKERRQ(ierr);
1672e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1673e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%d",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,value[0]);CHKERRQ(ierr);
1674e2446a98SMatthew Knepley     for (i=1; i<*n; i++) {
1675e55864a3SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%d",value[i]);CHKERRQ(ierr);
1676e2446a98SMatthew Knepley     }
1677e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
1678e2446a98SMatthew Knepley   }
1679e2446a98SMatthew Knepley   PetscFunctionReturn(0);
1680e2446a98SMatthew Knepley }
1681e2446a98SMatthew Knepley 
16828cc676e6SMatthew G Knepley /*@C
1683d1da0b69SBarry Smith    PetscOptionsViewer - Gets a viewer appropriate for the type indicated by the user
16848cc676e6SMatthew G Knepley 
16858cc676e6SMatthew G Knepley    Logically Collective on the communicator passed in PetscOptionsBegin()
16868cc676e6SMatthew G Knepley 
16878cc676e6SMatthew G Knepley    Input Parameters:
16888cc676e6SMatthew G Knepley +  opt - option name
16898cc676e6SMatthew G Knepley .  text - short string that describes the option
16908cc676e6SMatthew G Knepley -  man - manual page with additional information on option
16918cc676e6SMatthew G Knepley 
16928cc676e6SMatthew G Knepley    Output Parameter:
16938cc676e6SMatthew G Knepley +  viewer - the viewer
16948cc676e6SMatthew G Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
16958cc676e6SMatthew G Knepley 
16968cc676e6SMatthew G Knepley    Level: beginner
16978cc676e6SMatthew G Knepley 
16988cc676e6SMatthew G Knepley    Concepts: options database^has int
16998cc676e6SMatthew G Knepley 
17008cc676e6SMatthew G Knepley    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
17018cc676e6SMatthew G Knepley 
17025a7113b9SPatrick Sanan    See PetscOptionsGetViewer() for the format of the supplied viewer and its options
17038cc676e6SMatthew G Knepley 
1704c5929fdfSBarry Smith .seealso: PetscOptionsGetViewer(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(NULL,),
17058cc676e6SMatthew G Knepley           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
17068cc676e6SMatthew G Knepley           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
17078cc676e6SMatthew G Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
17088cc676e6SMatthew G Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
17098cc676e6SMatthew G Knepley           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1710a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
17118cc676e6SMatthew G Knepley @*/
17124416b707SBarry Smith PetscErrorCode  PetscOptionsViewer_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscViewer *viewer,PetscViewerFormat *format,PetscBool  *set)
17138cc676e6SMatthew G Knepley {
17148cc676e6SMatthew G Knepley   PetscErrorCode  ierr;
17154416b707SBarry Smith   PetscOptionItem amsopt;
17168cc676e6SMatthew G Knepley 
17178cc676e6SMatthew G Knepley   PetscFunctionBegin;
17181a1499c8SBarry Smith   if (!PetscOptionsObject->count) {
17194416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
172064facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
17215b02f95dSBarry Smith     ierr = PetscStrdup("",(char**)&amsopt->data);CHKERRQ(ierr);
17228cc676e6SMatthew G Knepley   }
1723e55864a3SBarry Smith   ierr = PetscOptionsGetViewer(PetscOptionsObject->comm,PetscOptionsObject->prefix,opt,viewer,format,set);CHKERRQ(ierr);
1724e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1725e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%s>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,"",text,ManSection(man));CHKERRQ(ierr);
17268cc676e6SMatthew G Knepley   }
17278cc676e6SMatthew G Knepley   PetscFunctionReturn(0);
17288cc676e6SMatthew G Knepley }
17298cc676e6SMatthew G Knepley 
173053acd3b1SBarry Smith 
173153acd3b1SBarry Smith /*@C
1732b52f573bSBarry Smith      PetscOptionsHead - Puts a heading before listing any more published options. Used, for example,
173353acd3b1SBarry Smith             in KSPSetFromOptions_GMRES().
173453acd3b1SBarry Smith 
17353f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
173653acd3b1SBarry Smith 
173753acd3b1SBarry Smith    Input Parameter:
173853acd3b1SBarry Smith .   head - the heading text
173953acd3b1SBarry Smith 
174053acd3b1SBarry Smith 
174153acd3b1SBarry Smith    Level: intermediate
174253acd3b1SBarry Smith 
174353acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
174453acd3b1SBarry Smith 
1745b52f573bSBarry Smith           Can be followed by a call to PetscOptionsTail() in the same function.
174653acd3b1SBarry Smith 
174753acd3b1SBarry Smith    Concepts: options database^subheading
174853acd3b1SBarry Smith 
1749c5929fdfSBarry Smith .seealso: PetscOptionsGetInt(NULL,), PetscOptionsGetReal(),
1750acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
175153acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
175253acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1753acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1754a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
175553acd3b1SBarry Smith @*/
17564416b707SBarry Smith PetscErrorCode  PetscOptionsHead(PetscOptionItems *PetscOptionsObject,const char head[])
175753acd3b1SBarry Smith {
175853acd3b1SBarry Smith   PetscErrorCode ierr;
175953acd3b1SBarry Smith 
176053acd3b1SBarry Smith   PetscFunctionBegin;
1761e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1762e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  %s\n",head);CHKERRQ(ierr);
176353acd3b1SBarry Smith   }
176453acd3b1SBarry Smith   PetscFunctionReturn(0);
176553acd3b1SBarry Smith }
176653acd3b1SBarry Smith 
176753acd3b1SBarry Smith 
176853acd3b1SBarry Smith 
176953acd3b1SBarry Smith 
177053acd3b1SBarry Smith 
177153acd3b1SBarry Smith 
1772