xref: /petsc/src/sys/objects/aoptions.c (revision 0e4c290a6c0eb17dbbd82d77cc66a07df1a900c0)
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 
432d747510SLisandro Dalcin   ierr = PetscOptionsHasHelp(PetscOptionsObject->options,&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 
111d083f849SBarry Smith     Collective
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 
16995452b02SPatrick Sanan     Notes:
17095452b02SPatrick Sanan     this isn't really practical, it is just to demonstrate the principle
171aee2cecaSBarry Smith 
1727781c08eSBarry Smith     A carriage return indicates no change from the default; but this like -ksp_monitor <stdout>  the default is actually not stdout the default
1737781c08eSBarry Smith     is to do nothing so to get it to use stdout you need to type stdout. This is kind of bug?
1747781c08eSBarry Smith 
175aee2cecaSBarry Smith     Bugs:
1767781c08eSBarry Smith +    All processes must traverse through the exact same set of option queries due to the call to PetscScanString()
1773cc1e11dSBarry Smith .    Internal strings have arbitrary length and string copies are not checked that they fit into string space
178aee2cecaSBarry Smith -    Only works for PetscInt == int, PetscReal == double etc
179aee2cecaSBarry Smith 
18095452b02SPatrick Sanan     Developer Notes:
18195452b02SPatrick Sanan     Normally the GUI that presents the options the user and retrieves the values would be running in a different
1823cc1e11dSBarry Smith      address space and communicating with the PETSc program
1833cc1e11dSBarry Smith 
184aee2cecaSBarry Smith */
1854416b707SBarry Smith PetscErrorCode PetscOptionsGetFromTextInput(PetscOptionItems *PetscOptionsObject)
1866356e834SBarry Smith {
1876356e834SBarry Smith   PetscErrorCode  ierr;
1884416b707SBarry Smith   PetscOptionItem next = PetscOptionsObject->next;
1896356e834SBarry Smith   char            str[512];
1907781c08eSBarry Smith   PetscBool       bid;
191a4404d99SBarry Smith   PetscReal       ir,*valr;
192330cf3c9SBarry Smith   PetscInt        *vald;
193330cf3c9SBarry Smith   size_t          i;
1946356e834SBarry Smith 
1952a409bb0SBarry Smith   PetscFunctionBegin;
196e55864a3SBarry Smith   ierr = (*PetscPrintf)(PETSC_COMM_WORLD,"%s -------------------------------------------------\n",PetscOptionsObject->title);CHKERRQ(ierr);
1976356e834SBarry Smith   while (next) {
1986356e834SBarry Smith     switch (next->type) {
1996356e834SBarry Smith     case OPTION_HEAD:
2006356e834SBarry Smith       break;
201e26ddf31SBarry Smith     case OPTION_INT_ARRAY:
202e55864a3SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",next->option+1);CHKERRQ(ierr);
203e26ddf31SBarry Smith       vald = (PetscInt*) next->data;
204e26ddf31SBarry Smith       for (i=0; i<next->arraylength; i++) {
205e26ddf31SBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"%d",vald[i]);CHKERRQ(ierr);
206e26ddf31SBarry Smith         if (i < next->arraylength-1) {
207e26ddf31SBarry Smith           ierr = PetscPrintf(PETSC_COMM_WORLD,",");CHKERRQ(ierr);
208e26ddf31SBarry Smith         }
209e26ddf31SBarry Smith       }
210e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,">: %s (%s) ",next->text,next->man);CHKERRQ(ierr);
211e26ddf31SBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
212e26ddf31SBarry Smith       if (str[0]) {
213e26ddf31SBarry Smith         PetscToken token;
214e26ddf31SBarry Smith         PetscInt   n=0,nmax = next->arraylength,*dvalue = (PetscInt*)next->data,start,end;
215e26ddf31SBarry Smith         size_t     len;
216e26ddf31SBarry Smith         char       *value;
217ace3abfcSBarry Smith         PetscBool  foundrange;
218e26ddf31SBarry Smith 
219e26ddf31SBarry Smith         next->set = PETSC_TRUE;
220e26ddf31SBarry Smith         value     = str;
221e26ddf31SBarry Smith         ierr      = PetscTokenCreate(value,',',&token);CHKERRQ(ierr);
222e26ddf31SBarry Smith         ierr      = PetscTokenFind(token,&value);CHKERRQ(ierr);
223e26ddf31SBarry Smith         while (n < nmax) {
224e26ddf31SBarry Smith           if (!value) break;
225e26ddf31SBarry Smith 
226e26ddf31SBarry Smith           /* look for form  d-D where d and D are integers */
227e26ddf31SBarry Smith           foundrange = PETSC_FALSE;
228e26ddf31SBarry Smith           ierr       = PetscStrlen(value,&len);CHKERRQ(ierr);
229e26ddf31SBarry Smith           if (value[0] == '-') i=2;
230e26ddf31SBarry Smith           else i=1;
231330cf3c9SBarry Smith           for (;i<len; i++) {
232e26ddf31SBarry Smith             if (value[i] == '-') {
233e32f2f54SBarry Smith               if (i == len-1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Error in %D-th array entry %s\n",n,value);
234e26ddf31SBarry Smith               value[i] = 0;
235cfbddea1SSatish Balay               ierr     = PetscOptionsStringToInt(value,&start);CHKERRQ(ierr);
236cfbddea1SSatish Balay               ierr     = PetscOptionsStringToInt(value+i+1,&end);CHKERRQ(ierr);
237e32f2f54SBarry 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);
238e32f2f54SBarry 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);
239e26ddf31SBarry Smith               for (; start<end; start++) {
240e26ddf31SBarry Smith                 *dvalue = start; dvalue++;n++;
241e26ddf31SBarry Smith               }
242e26ddf31SBarry Smith               foundrange = PETSC_TRUE;
243e26ddf31SBarry Smith               break;
244e26ddf31SBarry Smith             }
245e26ddf31SBarry Smith           }
246e26ddf31SBarry Smith           if (!foundrange) {
247cfbddea1SSatish Balay             ierr = PetscOptionsStringToInt(value,dvalue);CHKERRQ(ierr);
248e26ddf31SBarry Smith             dvalue++;
249e26ddf31SBarry Smith             n++;
250e26ddf31SBarry Smith           }
251e26ddf31SBarry Smith           ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
252e26ddf31SBarry Smith         }
2538c74ee41SBarry Smith         ierr = PetscTokenDestroy(&token);CHKERRQ(ierr);
254e26ddf31SBarry Smith       }
255e26ddf31SBarry Smith       break;
256e26ddf31SBarry Smith     case OPTION_REAL_ARRAY:
257e55864a3SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",next->option+1);CHKERRQ(ierr);
258e26ddf31SBarry Smith       valr = (PetscReal*) next->data;
259e26ddf31SBarry Smith       for (i=0; i<next->arraylength; i++) {
260e26ddf31SBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"%g",valr[i]);CHKERRQ(ierr);
261e26ddf31SBarry Smith         if (i < next->arraylength-1) {
262e26ddf31SBarry Smith           ierr = PetscPrintf(PETSC_COMM_WORLD,",");CHKERRQ(ierr);
263e26ddf31SBarry Smith         }
264e26ddf31SBarry Smith       }
265e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,">: %s (%s) ",next->text,next->man);CHKERRQ(ierr);
266e26ddf31SBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
267e26ddf31SBarry Smith       if (str[0]) {
268e26ddf31SBarry Smith         PetscToken token;
269e26ddf31SBarry Smith         PetscInt   n = 0,nmax = next->arraylength;
270e26ddf31SBarry Smith         PetscReal  *dvalue = (PetscReal*)next->data;
271e26ddf31SBarry Smith         char       *value;
272e26ddf31SBarry Smith 
273e26ddf31SBarry Smith         next->set = PETSC_TRUE;
274e26ddf31SBarry Smith         value     = str;
275e26ddf31SBarry Smith         ierr      = PetscTokenCreate(value,',',&token);CHKERRQ(ierr);
276e26ddf31SBarry Smith         ierr      = PetscTokenFind(token,&value);CHKERRQ(ierr);
277e26ddf31SBarry Smith         while (n < nmax) {
278e26ddf31SBarry Smith           if (!value) break;
279cfbddea1SSatish Balay           ierr = PetscOptionsStringToReal(value,dvalue);CHKERRQ(ierr);
280e26ddf31SBarry Smith           dvalue++;
281e26ddf31SBarry Smith           n++;
282e26ddf31SBarry Smith           ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
283e26ddf31SBarry Smith         }
2848c74ee41SBarry Smith         ierr = PetscTokenDestroy(&token);CHKERRQ(ierr);
285e26ddf31SBarry Smith       }
286e26ddf31SBarry Smith       break;
2876356e834SBarry Smith     case OPTION_INT:
288e55864a3SBarry 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);
2893fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
2903fc1eb6aSBarry Smith       if (str[0]) {
291d25d7f95SJed Brown #if defined(PETSC_SIZEOF_LONG_LONG)
292d25d7f95SJed Brown         long long lid;
293d25d7f95SJed Brown         sscanf(str,"%lld",&lid);
2941a1499c8SBarry 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);
295c272547aSJed Brown #else
296d25d7f95SJed Brown         long  lid;
297d25d7f95SJed Brown         sscanf(str,"%ld",&lid);
2981a1499c8SBarry 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);
299c272547aSJed Brown #endif
300a297a907SKarl Rupp 
301d25d7f95SJed Brown         next->set = PETSC_TRUE;
302d25d7f95SJed Brown         *((PetscInt*)next->data) = (PetscInt)lid;
303aee2cecaSBarry Smith       }
304aee2cecaSBarry Smith       break;
305aee2cecaSBarry Smith     case OPTION_REAL:
306e55864a3SBarry 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);
3073fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3083fc1eb6aSBarry Smith       if (str[0]) {
309ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
310a4404d99SBarry Smith         sscanf(str,"%e",&ir);
311570b7f6dSBarry Smith #elif defined(PETSC_USE_REAL___FP16)
312570b7f6dSBarry Smith         float irtemp;
313570b7f6dSBarry Smith         sscanf(str,"%e",&irtemp);
314570b7f6dSBarry Smith         ir = irtemp;
315ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL_DOUBLE)
316aee2cecaSBarry Smith         sscanf(str,"%le",&ir);
317ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL___FLOAT128)
318d9822059SBarry Smith         ir = strtoflt128(str,0);
319d9822059SBarry Smith #else
320513dbe71SLisandro Dalcin         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Unknown scalar type");
321a4404d99SBarry Smith #endif
322aee2cecaSBarry Smith         next->set                 = PETSC_TRUE;
323aee2cecaSBarry Smith         *((PetscReal*)next->data) = ir;
324aee2cecaSBarry Smith       }
325aee2cecaSBarry Smith       break;
3267781c08eSBarry Smith     case OPTION_BOOL:
32783355fc5SBarry 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);
3287781c08eSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3297781c08eSBarry Smith       if (str[0]) {
3307781c08eSBarry Smith         ierr = PetscOptionsStringToBool(str,&bid);CHKERRQ(ierr);
3317781c08eSBarry Smith         next->set = PETSC_TRUE;
3327781c08eSBarry Smith         *((PetscBool*)next->data) = bid;
3337781c08eSBarry Smith       }
3347781c08eSBarry Smith       break;
335aee2cecaSBarry Smith     case OPTION_STRING:
336e55864a3SBarry 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);
3373fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3383fc1eb6aSBarry Smith       if (str[0]) {
339aee2cecaSBarry Smith         next->set = PETSC_TRUE;
34064facd6cSBarry Smith         /* must use system malloc since SAWs may free this */
3415b02f95dSBarry Smith         ierr = PetscStrdup(str,(char**)&next->data);CHKERRQ(ierr);
3426356e834SBarry Smith       }
3436356e834SBarry Smith       break;
344a264d7a6SBarry Smith     case OPTION_FLIST:
34544ef3d73SBarry Smith       ierr = PetscFunctionListPrintTypes(PETSC_COMM_WORLD,stdout,PetscOptionsObject->prefix,next->option,next->text,next->man,next->flist,(char*)next->data,(char*)next->data);CHKERRQ(ierr);
3463cc1e11dSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3473cc1e11dSBarry Smith       if (str[0]) {
348e55864a3SBarry Smith         PetscOptionsObject->changedmethod = PETSC_TRUE;
3493cc1e11dSBarry Smith         next->set = PETSC_TRUE;
35064facd6cSBarry Smith         /* must use system malloc since SAWs may free this */
3515b02f95dSBarry Smith         ierr = PetscStrdup(str,(char**)&next->data);CHKERRQ(ierr);
3523cc1e11dSBarry Smith       }
3533cc1e11dSBarry Smith       break;
354b432afdaSMatthew Knepley     default:
355b432afdaSMatthew Knepley       break;
3566356e834SBarry Smith     }
3576356e834SBarry Smith     next = next->next;
3586356e834SBarry Smith   }
3596356e834SBarry Smith   PetscFunctionReturn(0);
3606356e834SBarry Smith }
3616356e834SBarry Smith 
362e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
363e04113cfSBarry Smith #include <petscviewersaws.h>
364d5649816SBarry Smith 
365d5649816SBarry Smith static int count = 0;
366d5649816SBarry Smith 
367e04113cfSBarry Smith PetscErrorCode PetscOptionsSAWsDestroy(void)
368d5649816SBarry Smith {
3692657e9d9SBarry Smith   PetscFunctionBegin;
370d5649816SBarry Smith   PetscFunctionReturn(0);
371d5649816SBarry Smith }
372d5649816SBarry Smith 
3739c1e0ce8SBarry Smith static const char *OptionsHeader = "<head>\n"
374a8d69d7bSBarry Smith                                    "<script type=\"text/javascript\" src=\"https://www.mcs.anl.gov/research/projects/saws/js/jquery-1.9.1.js\"></script>\n"
375a8d69d7bSBarry Smith                                    "<script type=\"text/javascript\" src=\"https://www.mcs.anl.gov/research/projects/saws/js/SAWs.js\"></script>\n"
376d1fc0251SBarry Smith                                    "<script type=\"text/javascript\" src=\"js/PETSc.js\"></script>\n"
37764bbc9efSBarry Smith                                    "<script>\n"
37864bbc9efSBarry Smith                                       "jQuery(document).ready(function() {\n"
37964bbc9efSBarry Smith                                          "PETSc.getAndDisplayDirectory(null,\"#variablesInfo\")\n"
38064bbc9efSBarry Smith                                       "})\n"
38164bbc9efSBarry Smith                                   "</script>\n"
38264bbc9efSBarry Smith                                   "</head>\n";
3831423471aSBarry Smith 
3841423471aSBarry Smith /*  Determines the size and style of the scroll region where PETSc options selectable from users are displayed */
3851423471aSBarry Smith static const char *OptionsBodyBottom = "<div id=\"variablesInfo\" style=\"background-color:lightblue;height:auto;max-height:500px;overflow:scroll;\"></div>\n<br>\n</body>";
38664bbc9efSBarry Smith 
387b3506946SBarry Smith /*
3887781c08eSBarry Smith     PetscOptionsSAWsInput - Presents all the PETSc Options processed by the program so the user may change them at runtime using the SAWs
389b3506946SBarry Smith 
390b3506946SBarry Smith     Bugs:
391b3506946SBarry Smith +    All processes must traverse through the exact same set of option queries do to the call to PetscScanString()
392b3506946SBarry Smith .    Internal strings have arbitrary length and string copies are not checked that they fit into string space
393b3506946SBarry Smith -    Only works for PetscInt == int, PetscReal == double etc
394b3506946SBarry Smith 
395b3506946SBarry Smith 
396b3506946SBarry Smith */
3974416b707SBarry Smith PetscErrorCode PetscOptionsSAWsInput(PetscOptionItems *PetscOptionsObject)
398b3506946SBarry Smith {
399b3506946SBarry Smith   PetscErrorCode  ierr;
4004416b707SBarry Smith   PetscOptionItem next     = PetscOptionsObject->next;
401d5649816SBarry Smith   static int      mancount = 0;
402b3506946SBarry Smith   char            options[16];
4037aab2a10SBarry Smith   PetscBool       changedmethod = PETSC_FALSE;
404a23eb982SSurtai Han   PetscBool       stopasking    = PETSC_FALSE;
40588a9752cSBarry Smith   char            manname[16],textname[16];
4062657e9d9SBarry Smith   char            dir[1024];
407b3506946SBarry Smith 
4082a409bb0SBarry Smith   PetscFunctionBegin;
409b3506946SBarry Smith   /* the next line is a bug, this will only work if all processors are here, the comm passed in is ignored!!! */
410b3506946SBarry Smith   sprintf(options,"Options_%d",count++);
411a297a907SKarl Rupp 
4127325c4a4SBarry Smith   PetscOptionsObject->pprefix = PetscOptionsObject->prefix; /* SAWs will change this, so cannot pass prefix directly */
4131bc75a8dSBarry Smith 
414d50831c4SBarry Smith   ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s","_title");CHKERRQ(ierr);
41583355fc5SBarry Smith   PetscStackCallSAWs(SAWs_Register,(dir,&PetscOptionsObject->title,1,SAWs_READ,SAWs_STRING));
4167781c08eSBarry Smith   ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s","prefix");CHKERRQ(ierr);
41783355fc5SBarry Smith   PetscStackCallSAWs(SAWs_Register,(dir,&PetscOptionsObject->pprefix,1,SAWs_READ,SAWs_STRING));
4182657e9d9SBarry Smith   PetscStackCallSAWs(SAWs_Register,("/PETSc/Options/ChangedMethod",&changedmethod,1,SAWs_WRITE,SAWs_BOOLEAN));
419a23eb982SSurtai Han   PetscStackCallSAWs(SAWs_Register,("/PETSc/Options/StopAsking",&stopasking,1,SAWs_WRITE,SAWs_BOOLEAN));
420b3506946SBarry Smith 
421b3506946SBarry Smith   while (next) {
422d50831c4SBarry Smith     sprintf(manname,"_man_%d",mancount);
4232657e9d9SBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",manname);CHKERRQ(ierr);
4247781c08eSBarry Smith     PetscStackCallSAWs(SAWs_Register,(dir,&next->man,1,SAWs_READ,SAWs_STRING));
425d50831c4SBarry Smith     sprintf(textname,"_text_%d",mancount++);
4267781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",textname);CHKERRQ(ierr);
4277781c08eSBarry Smith     PetscStackCallSAWs(SAWs_Register,(dir,&next->text,1,SAWs_READ,SAWs_STRING));
4289f32e415SBarry Smith 
429b3506946SBarry Smith     switch (next->type) {
430b3506946SBarry Smith     case OPTION_HEAD:
431b3506946SBarry Smith       break;
432b3506946SBarry Smith     case OPTION_INT_ARRAY:
4337781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4342657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_INT));
435b3506946SBarry Smith       break;
436b3506946SBarry Smith     case OPTION_REAL_ARRAY:
4377781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4382657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_DOUBLE));
439b3506946SBarry Smith       break;
440b3506946SBarry Smith     case OPTION_INT:
4417781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4422657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_INT));
443b3506946SBarry Smith       break;
444b3506946SBarry Smith     case OPTION_REAL:
4457781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4462657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_DOUBLE));
447b3506946SBarry Smith       break;
4487781c08eSBarry Smith     case OPTION_BOOL:
4497781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4502657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_BOOLEAN));
4511ae3d29cSBarry Smith       break;
4527781c08eSBarry Smith     case OPTION_BOOL_ARRAY:
4537781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4542657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_BOOLEAN));
45571f08665SBarry Smith       break;
456b3506946SBarry Smith     case OPTION_STRING:
4577781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4587781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
4591ae3d29cSBarry Smith       break;
4601ae3d29cSBarry Smith     case OPTION_STRING_ARRAY:
4617781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4622657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_STRING));
463b3506946SBarry Smith       break;
464a264d7a6SBarry Smith     case OPTION_FLIST:
465a264d7a6SBarry Smith       {
466a264d7a6SBarry Smith       PetscInt ntext;
4677781c08eSBarry Smith       ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4687781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
469a264d7a6SBarry Smith       ierr = PetscFunctionListGet(next->flist,(const char***)&next->edata,&ntext);CHKERRQ(ierr);
470a264d7a6SBarry Smith       PetscStackCallSAWs(SAWs_Set_Legal_Variable_Values,(dir,ntext,next->edata));
471a264d7a6SBarry Smith       }
472a264d7a6SBarry Smith       break;
4731ae3d29cSBarry Smith     case OPTION_ELIST:
474a264d7a6SBarry Smith       {
475a264d7a6SBarry Smith       PetscInt ntext = next->nlist;
4767781c08eSBarry Smith       ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4777781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
478ead66b60SBarry Smith       ierr = PetscMalloc1((ntext+1),(char***)&next->edata);CHKERRQ(ierr);
4791ae3d29cSBarry Smith       ierr = PetscMemcpy(next->edata,next->list,ntext*sizeof(char*));CHKERRQ(ierr);
480a264d7a6SBarry Smith       PetscStackCallSAWs(SAWs_Set_Legal_Variable_Values,(dir,ntext,next->edata));
481a264d7a6SBarry Smith       }
482a264d7a6SBarry Smith       break;
483b3506946SBarry Smith     default:
484b3506946SBarry Smith       break;
485b3506946SBarry Smith     }
486b3506946SBarry Smith     next = next->next;
487b3506946SBarry Smith   }
488b3506946SBarry Smith 
489b3506946SBarry Smith   /* wait until accessor has unlocked the memory */
49064bbc9efSBarry Smith   PetscStackCallSAWs(SAWs_Push_Header,("index.html",OptionsHeader));
49164bbc9efSBarry Smith   PetscStackCallSAWs(SAWs_Push_Body,("index.html",2,OptionsBodyBottom));
4927aab2a10SBarry Smith   ierr = PetscSAWsBlock();CHKERRQ(ierr);
49364bbc9efSBarry Smith   PetscStackCallSAWs(SAWs_Pop_Header,("index.html"));
49464bbc9efSBarry Smith   PetscStackCallSAWs(SAWs_Pop_Body,("index.html",2));
495b3506946SBarry Smith 
49688a9752cSBarry Smith   /* determine if any values have been set in GUI */
49783355fc5SBarry Smith   next = PetscOptionsObject->next;
49888a9752cSBarry Smith   while (next) {
49988a9752cSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
500f7b25cf6SBarry Smith     PetscStackCallSAWs(SAWs_Selected,(dir,(int*)&next->set));
50188a9752cSBarry Smith     next = next->next;
50288a9752cSBarry Smith   }
50388a9752cSBarry Smith 
504b3506946SBarry Smith   /* reset counter to -2; this updates the screen with the new options for the selected method */
505f7b25cf6SBarry Smith   if (changedmethod) PetscOptionsObject->count = -2;
506b3506946SBarry Smith 
507a23eb982SSurtai Han   if (stopasking) {
508a23eb982SSurtai Han     PetscOptionsPublish      = PETSC_FALSE;
50912655325SBarry Smith     PetscOptionsObject->count = 0;//do not ask for same thing again
510a23eb982SSurtai Han   }
511a23eb982SSurtai Han 
5129a492a5cSBarry Smith   PetscStackCallSAWs(SAWs_Delete,("/PETSc/Options"));
513b3506946SBarry Smith   PetscFunctionReturn(0);
514b3506946SBarry Smith }
515b3506946SBarry Smith #endif
516b3506946SBarry Smith 
5174416b707SBarry Smith PetscErrorCode PetscOptionsEnd_Private(PetscOptionItems *PetscOptionsObject)
51853acd3b1SBarry Smith {
51953acd3b1SBarry Smith   PetscErrorCode  ierr;
5204416b707SBarry Smith   PetscOptionItem last;
5216356e834SBarry Smith   char            option[256],value[1024],tmp[32];
522330cf3c9SBarry Smith   size_t          j;
52353acd3b1SBarry Smith 
52453acd3b1SBarry Smith   PetscFunctionBegin;
52583355fc5SBarry Smith   if (PetscOptionsObject->next) {
52683355fc5SBarry Smith     if (!PetscOptionsObject->count) {
527a264d7a6SBarry Smith #if defined(PETSC_HAVE_SAWS)
528f7b25cf6SBarry Smith       ierr = PetscOptionsSAWsInput(PetscOptionsObject);CHKERRQ(ierr);
529b3506946SBarry Smith #else
530e55864a3SBarry Smith       ierr = PetscOptionsGetFromTextInput(PetscOptionsObject);CHKERRQ(ierr);
531b3506946SBarry Smith #endif
532aee2cecaSBarry Smith     }
533aee2cecaSBarry Smith   }
5346356e834SBarry Smith 
535e55864a3SBarry Smith   ierr = PetscFree(PetscOptionsObject->title);CHKERRQ(ierr);
5366356e834SBarry Smith 
537e26ddf31SBarry Smith   /* reset counter to -2; this updates the screen with the new options for the selected method */
538e55864a3SBarry Smith   if (PetscOptionsObject->changedmethod) PetscOptionsObject->count = -2;
5397a72a596SBarry Smith   /* reset alreadyprinted flag */
540e55864a3SBarry Smith   PetscOptionsObject->alreadyprinted = PETSC_FALSE;
541e55864a3SBarry Smith   if (PetscOptionsObject->object) PetscOptionsObject->object->optionsprinted = PETSC_TRUE;
542e55864a3SBarry Smith   PetscOptionsObject->object = NULL;
54353acd3b1SBarry Smith 
544e55864a3SBarry Smith   while (PetscOptionsObject->next) {
545e55864a3SBarry Smith     if (PetscOptionsObject->next->set) {
546e55864a3SBarry Smith       if (PetscOptionsObject->prefix) {
54753acd3b1SBarry Smith         ierr = PetscStrcpy(option,"-");CHKERRQ(ierr);
548e55864a3SBarry Smith         ierr = PetscStrcat(option,PetscOptionsObject->prefix);CHKERRQ(ierr);
549e55864a3SBarry Smith         ierr = PetscStrcat(option,PetscOptionsObject->next->option+1);CHKERRQ(ierr);
5506356e834SBarry Smith       } else {
551e55864a3SBarry Smith         ierr = PetscStrcpy(option,PetscOptionsObject->next->option);CHKERRQ(ierr);
5526356e834SBarry Smith       }
5536356e834SBarry Smith 
554e55864a3SBarry Smith       switch (PetscOptionsObject->next->type) {
5556356e834SBarry Smith       case OPTION_HEAD:
5566356e834SBarry Smith         break;
5576356e834SBarry Smith       case OPTION_INT_ARRAY:
558e55864a3SBarry Smith         sprintf(value,"%d",(int)((PetscInt*)PetscOptionsObject->next->data)[0]);
559e55864a3SBarry Smith         for (j=1; j<PetscOptionsObject->next->arraylength; j++) {
560e55864a3SBarry Smith           sprintf(tmp,"%d",(int)((PetscInt*)PetscOptionsObject->next->data)[j]);
5616356e834SBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5626356e834SBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5636356e834SBarry Smith         }
5646356e834SBarry Smith         break;
5656356e834SBarry Smith       case OPTION_INT:
566e55864a3SBarry Smith         sprintf(value,"%d",(int) *(PetscInt*)PetscOptionsObject->next->data);
5676356e834SBarry Smith         break;
5686356e834SBarry Smith       case OPTION_REAL:
569e55864a3SBarry Smith         sprintf(value,"%g",(double) *(PetscReal*)PetscOptionsObject->next->data);
5706356e834SBarry Smith         break;
5716356e834SBarry Smith       case OPTION_REAL_ARRAY:
572e55864a3SBarry Smith         sprintf(value,"%g",(double)((PetscReal*)PetscOptionsObject->next->data)[0]);
573e55864a3SBarry Smith         for (j=1; j<PetscOptionsObject->next->arraylength; j++) {
574e55864a3SBarry Smith           sprintf(tmp,"%g",(double)((PetscReal*)PetscOptionsObject->next->data)[j]);
5756356e834SBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5766356e834SBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5776356e834SBarry Smith         }
5786356e834SBarry Smith         break;
579050cccc3SHong Zhang       case OPTION_SCALAR_ARRAY:
58095f3a755SJose E. Roman         sprintf(value,"%g+%gi",(double)PetscRealPart(((PetscScalar*)PetscOptionsObject->next->data)[0]),(double)PetscImaginaryPart(((PetscScalar*)PetscOptionsObject->next->data)[0]));
581050cccc3SHong Zhang         for (j=1; j<PetscOptionsObject->next->arraylength; j++) {
58295f3a755SJose E. Roman           sprintf(tmp,"%g+%gi",(double)PetscRealPart(((PetscScalar*)PetscOptionsObject->next->data)[j]),(double)PetscImaginaryPart(((PetscScalar*)PetscOptionsObject->next->data)[j]));
583050cccc3SHong Zhang           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
584050cccc3SHong Zhang           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
585050cccc3SHong Zhang         }
586050cccc3SHong Zhang         break;
5877781c08eSBarry Smith       case OPTION_BOOL:
588e55864a3SBarry Smith         sprintf(value,"%d",*(int*)PetscOptionsObject->next->data);
5896356e834SBarry Smith         break;
5907781c08eSBarry Smith       case OPTION_BOOL_ARRAY:
591e55864a3SBarry Smith         sprintf(value,"%d",(int)((PetscBool*)PetscOptionsObject->next->data)[0]);
592e55864a3SBarry Smith         for (j=1; j<PetscOptionsObject->next->arraylength; j++) {
593e55864a3SBarry Smith           sprintf(tmp,"%d",(int)((PetscBool*)PetscOptionsObject->next->data)[j]);
5941ae3d29cSBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5951ae3d29cSBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5961ae3d29cSBarry Smith         }
5971ae3d29cSBarry Smith         break;
598a264d7a6SBarry Smith       case OPTION_FLIST:
5996991f827SBarry Smith         ierr = PetscStrcpy(value,(char*)PetscOptionsObject->next->data);CHKERRQ(ierr);
6006991f827SBarry Smith         break;
6011ae3d29cSBarry Smith       case OPTION_ELIST:
602e55864a3SBarry Smith         ierr = PetscStrcpy(value,(char*)PetscOptionsObject->next->data);CHKERRQ(ierr);
6036356e834SBarry Smith         break;
6041ae3d29cSBarry Smith       case OPTION_STRING:
605e55864a3SBarry Smith         ierr = PetscStrcpy(value,(char*)PetscOptionsObject->next->data);CHKERRQ(ierr);
606d51da6bfSBarry Smith         break;
6071ae3d29cSBarry Smith       case OPTION_STRING_ARRAY:
608e55864a3SBarry Smith         sprintf(value,"%s",((char**)PetscOptionsObject->next->data)[0]);
609e55864a3SBarry Smith         for (j=1; j<PetscOptionsObject->next->arraylength; j++) {
610e55864a3SBarry Smith           sprintf(tmp,"%s",((char**)PetscOptionsObject->next->data)[j]);
6111ae3d29cSBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
6121ae3d29cSBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
6131ae3d29cSBarry Smith         }
6146356e834SBarry Smith         break;
6156356e834SBarry Smith       }
616c5929fdfSBarry Smith       ierr = PetscOptionsSetValue(PetscOptionsObject->options,option,value);CHKERRQ(ierr);
6176356e834SBarry Smith     }
6186991f827SBarry Smith     if (PetscOptionsObject->next->type == OPTION_ELIST) {
6196991f827SBarry Smith       ierr = PetscStrNArrayDestroy(PetscOptionsObject->next->nlist,(char ***)&PetscOptionsObject->next->list);CHKERRQ(ierr);
6206991f827SBarry Smith     }
621e55864a3SBarry Smith     ierr   = PetscFree(PetscOptionsObject->next->text);CHKERRQ(ierr);
622e55864a3SBarry Smith     ierr   = PetscFree(PetscOptionsObject->next->option);CHKERRQ(ierr);
623e55864a3SBarry Smith     ierr   = PetscFree(PetscOptionsObject->next->man);CHKERRQ(ierr);
624e55864a3SBarry Smith     ierr   = PetscFree(PetscOptionsObject->next->edata);CHKERRQ(ierr);
625c979a496SBarry Smith 
62683355fc5SBarry Smith     if ((PetscOptionsObject->next->type == OPTION_STRING) || (PetscOptionsObject->next->type == OPTION_FLIST) || (PetscOptionsObject->next->type == OPTION_ELIST)){
62783355fc5SBarry Smith       free(PetscOptionsObject->next->data);
628c979a496SBarry Smith     } else {
62983355fc5SBarry Smith       ierr   = PetscFree(PetscOptionsObject->next->data);CHKERRQ(ierr);
630c979a496SBarry Smith     }
6317781c08eSBarry Smith 
63283355fc5SBarry Smith     last                    = PetscOptionsObject->next;
63383355fc5SBarry Smith     PetscOptionsObject->next = PetscOptionsObject->next->next;
6346356e834SBarry Smith     ierr                    = PetscFree(last);CHKERRQ(ierr);
6356356e834SBarry Smith   }
636f59f755dSBarry Smith   ierr = PetscFree(PetscOptionsObject->prefix);CHKERRQ(ierr);
637e55864a3SBarry Smith   PetscOptionsObject->next = 0;
63853acd3b1SBarry Smith   PetscFunctionReturn(0);
63953acd3b1SBarry Smith }
64053acd3b1SBarry Smith 
64153acd3b1SBarry Smith /*@C
64253acd3b1SBarry Smith    PetscOptionsEnum - Gets the enum value for a particular option in the database.
64353acd3b1SBarry Smith 
6443f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
64553acd3b1SBarry Smith 
64653acd3b1SBarry Smith    Input Parameters:
64753acd3b1SBarry Smith +  opt - option name
64853acd3b1SBarry Smith .  text - short string that describes the option
64953acd3b1SBarry Smith .  man - manual page with additional information on option
65053acd3b1SBarry Smith .  list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null
6510fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
6520fdccdaeSBarry Smith $                 PetscOptionsEnum(..., obj->value,&object->value,...) or
6530fdccdaeSBarry Smith $                 value = defaultvalue
6540fdccdaeSBarry Smith $                 PetscOptionsEnum(..., value,&value,&flg);
6550fdccdaeSBarry Smith $                 if (flg) {
65653acd3b1SBarry Smith 
65753acd3b1SBarry Smith    Output Parameter:
65853acd3b1SBarry Smith +  value - the  value to return
659b32e0204SMatthew G Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
66053acd3b1SBarry Smith 
66153acd3b1SBarry Smith    Level: beginner
66253acd3b1SBarry Smith 
66395452b02SPatrick Sanan    Notes:
66495452b02SPatrick Sanan     Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
66553acd3b1SBarry Smith 
66653acd3b1SBarry Smith           list is usually something like PCASMTypes or some other predefined list of enum names
66753acd3b1SBarry Smith 
6682efd9cb1SBarry Smith           If the user does not supply the option at all value is NOT changed. Thus
6692efd9cb1SBarry Smith           you should ALWAYS initialize value if you access it without first checking if the set flag is true.
6702efd9cb1SBarry Smith 
671989712b9SBarry Smith           The default/currentvalue passed into this routine does not get transferred to the output value variable automatically.
672989712b9SBarry Smith 
67353acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
674acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
675acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
67653acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
67753acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
678acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
679a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
68053acd3b1SBarry Smith @*/
6814416b707SBarry 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)
68253acd3b1SBarry Smith {
68353acd3b1SBarry Smith   PetscErrorCode ierr;
68453acd3b1SBarry Smith   PetscInt       ntext = 0;
685aa5bb8c0SSatish Balay   PetscInt       tval;
686ace3abfcSBarry Smith   PetscBool      tflg;
68753acd3b1SBarry Smith 
68853acd3b1SBarry Smith   PetscFunctionBegin;
68953acd3b1SBarry Smith   while (list[ntext++]) {
690e32f2f54SBarry Smith     if (ntext > 50) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument appears to be wrong or have more than 50 entries");
69153acd3b1SBarry Smith   }
692e32f2f54SBarry Smith   if (ntext < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument must have at least two entries: typename and type prefix");
69353acd3b1SBarry Smith   ntext -= 3;
694e55864a3SBarry Smith   ierr   = PetscOptionsEList_Private(PetscOptionsObject,opt,text,man,list,ntext,list[currentvalue],&tval,&tflg);CHKERRQ(ierr);
695aa5bb8c0SSatish Balay   /* with PETSC_USE_64BIT_INDICES sizeof(PetscInt) != sizeof(PetscEnum) */
696aa5bb8c0SSatish Balay   if (tflg) *value = (PetscEnum)tval;
697aa5bb8c0SSatish Balay   if (set)  *set   = tflg;
69853acd3b1SBarry Smith   PetscFunctionReturn(0);
69953acd3b1SBarry Smith }
70053acd3b1SBarry Smith 
701d3e47460SLisandro Dalcin /*@C
702d3e47460SLisandro Dalcin    PetscOptionsEnumArray - Gets an array of enum values for a particular
703d3e47460SLisandro Dalcin    option in the database.
704d3e47460SLisandro Dalcin 
705d3e47460SLisandro Dalcin    Logically Collective on the communicator passed in PetscOptionsBegin()
706d3e47460SLisandro Dalcin 
707d3e47460SLisandro Dalcin    Input Parameters:
708d3e47460SLisandro Dalcin +  opt - the option one is seeking
709d3e47460SLisandro Dalcin .  text - short string describing option
710d3e47460SLisandro Dalcin .  man - manual page for option
71122d58ec6SMatthew G. Knepley .  list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null
712d3e47460SLisandro Dalcin -  n - maximum number of values
713d3e47460SLisandro Dalcin 
714d3e47460SLisandro Dalcin    Output Parameter:
715d3e47460SLisandro Dalcin +  value - location to copy values
716d3e47460SLisandro Dalcin .  n - actual number of values found
717d3e47460SLisandro Dalcin -  set - PETSC_TRUE if found, else PETSC_FALSE
718d3e47460SLisandro Dalcin 
719d3e47460SLisandro Dalcin    Level: beginner
720d3e47460SLisandro Dalcin 
721d3e47460SLisandro Dalcin    Notes:
722d3e47460SLisandro Dalcin    The array must be passed as a comma separated list.
723d3e47460SLisandro Dalcin 
724d3e47460SLisandro Dalcin    There must be no intervening spaces between the values.
725d3e47460SLisandro Dalcin 
726d3e47460SLisandro Dalcin    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
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 
78695452b02SPatrick Sanan    Notes:
78795452b02SPatrick Sanan     If the user does not supply the option at all value is NOT changed. Thus
7882efd9cb1SBarry Smith           you should ALWAYS initialize value if you access it without first checking if the set flag is true.
7892efd9cb1SBarry Smith 
790989712b9SBarry Smith           The default/currentvalue passed into this routine does not get transferred to the output value variable automatically.
791989712b9SBarry Smith 
79253acd3b1SBarry Smith    Level: beginner
79353acd3b1SBarry Smith 
79495452b02SPatrick Sanan    Notes:
79595452b02SPatrick Sanan     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   }
82244ef3d73SBarry Smith   ierr = PetscOptionsGetInt(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,&wasset);CHKERRQ(ierr);
82344ef3d73SBarry Smith   if (set) *set = wasset;
824e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
825*0e4c290aSBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <now %D : formerly %D>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,wasset && value ? *value : currentvalue,currentvalue,text,ManSection(man));CHKERRQ(ierr);
82653acd3b1SBarry Smith   }
82753acd3b1SBarry Smith   PetscFunctionReturn(0);
82853acd3b1SBarry Smith }
82953acd3b1SBarry Smith 
83053acd3b1SBarry Smith /*@C
83153acd3b1SBarry Smith    PetscOptionsString - Gets the string value for a particular option in the database.
83253acd3b1SBarry Smith 
8333f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
83453acd3b1SBarry Smith 
83553acd3b1SBarry Smith    Input Parameters:
83653acd3b1SBarry Smith +  opt - option name
83753acd3b1SBarry Smith .  text - short string that describes the option
83853acd3b1SBarry Smith .  man - manual page with additional information on option
8390fdccdaeSBarry Smith .  currentvalue - the current value; caller is responsible for setting this value correctly. This is not used to set value
840bcbf2dc5SJed Brown -  len - length of the result string including null terminator
84153acd3b1SBarry Smith 
84253acd3b1SBarry Smith    Output Parameter:
84353acd3b1SBarry Smith +  value - the value to return
84453acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
84553acd3b1SBarry Smith 
84653acd3b1SBarry Smith    Level: beginner
84753acd3b1SBarry Smith 
84895452b02SPatrick Sanan    Notes:
84995452b02SPatrick Sanan     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 
856989712b9SBarry Smith           The default/currentvalue passed into this routine does not get transferred to the output value variable automatically.
857989712b9SBarry Smith 
8582efd9cb1SBarry Smith 
85989a13869SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
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;
87144ef3d73SBarry Smith   PetscBool       lset;
87253acd3b1SBarry Smith 
87353acd3b1SBarry Smith   PetscFunctionBegin;
8741a1499c8SBarry Smith   if (!PetscOptionsObject->count) {
8754416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
87664facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
8770fdccdaeSBarry Smith     ierr = PetscStrdup(currentvalue ? currentvalue : "",(char**)&amsopt->data);CHKERRQ(ierr);
878af6d86caSBarry Smith   }
87944ef3d73SBarry Smith   ierr = PetscOptionsGetString(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,len,&lset);CHKERRQ(ierr);
88044ef3d73SBarry Smith   if (set) *set = lset;
881e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
882*0e4c290aSBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <now %s : formerly %s>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,lset && value ? value : currentvalue,currentvalue,text,ManSection(man));CHKERRQ(ierr);
88353acd3b1SBarry Smith   }
88453acd3b1SBarry Smith   PetscFunctionReturn(0);
88553acd3b1SBarry Smith }
88653acd3b1SBarry Smith 
88753acd3b1SBarry Smith /*@C
88853acd3b1SBarry Smith    PetscOptionsReal - Gets the PetscReal value for a particular option in the database.
88953acd3b1SBarry Smith 
8903f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
89153acd3b1SBarry Smith 
89253acd3b1SBarry Smith    Input Parameters:
89353acd3b1SBarry Smith +  opt - option name
89453acd3b1SBarry Smith .  text - short string that describes the option
89553acd3b1SBarry Smith .  man - manual page with additional information on option
8960fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
8970fdccdaeSBarry Smith $                 PetscOptionsReal(..., obj->value,&object->value,...) or
8980fdccdaeSBarry Smith $                 value = defaultvalue
8990fdccdaeSBarry Smith $                 PetscOptionsReal(..., value,&value,&flg);
9000fdccdaeSBarry Smith $                 if (flg) {
90153acd3b1SBarry Smith 
90253acd3b1SBarry Smith    Output Parameter:
90353acd3b1SBarry Smith +  value - the value to return
90453acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
90553acd3b1SBarry Smith 
90695452b02SPatrick Sanan    Notes:
90795452b02SPatrick Sanan     If the user does not supply the option at all value is NOT changed. Thus
9082efd9cb1SBarry Smith           you should ALWAYS initialize value if you access it without first checking if the set flag is true.
9092efd9cb1SBarry Smith 
910989712b9SBarry Smith           The default/currentvalue passed into this routine does not get transferred to the output value variable automatically.
911989712b9SBarry Smith 
91253acd3b1SBarry Smith    Level: beginner
91353acd3b1SBarry Smith 
91495452b02SPatrick Sanan    Notes:
91595452b02SPatrick Sanan     Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
91653acd3b1SBarry Smith 
91789a13869SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
918acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
919acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
92053acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
92153acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
922acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
923a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
92453acd3b1SBarry Smith @*/
9254416b707SBarry Smith PetscErrorCode  PetscOptionsReal_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscReal currentvalue,PetscReal *value,PetscBool  *set)
92653acd3b1SBarry Smith {
92753acd3b1SBarry Smith   PetscErrorCode  ierr;
9284416b707SBarry Smith   PetscOptionItem amsopt;
92944ef3d73SBarry Smith   PetscBool       lset;
93053acd3b1SBarry Smith 
93153acd3b1SBarry Smith   PetscFunctionBegin;
932e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
9334416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_REAL,&amsopt);CHKERRQ(ierr);
934538aa990SBarry Smith     ierr = PetscMalloc(sizeof(PetscReal),&amsopt->data);CHKERRQ(ierr);
935a297a907SKarl Rupp 
9360fdccdaeSBarry Smith     *(PetscReal*)amsopt->data = currentvalue;
937538aa990SBarry Smith   }
93844ef3d73SBarry Smith   ierr = PetscOptionsGetReal(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,&lset);CHKERRQ(ierr);
93944ef3d73SBarry Smith   if (set) *set = lset;
9401a1499c8SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
941*0e4c290aSBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%g : %g>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,lset && value ? (double)*value : (double) currentvalue,(double)currentvalue,text,ManSection(man));CHKERRQ(ierr);
94253acd3b1SBarry Smith   }
94353acd3b1SBarry Smith   PetscFunctionReturn(0);
94453acd3b1SBarry Smith }
94553acd3b1SBarry Smith 
94653acd3b1SBarry Smith /*@C
94753acd3b1SBarry Smith    PetscOptionsScalar - Gets the scalar value for a particular option in the database.
94853acd3b1SBarry Smith 
9493f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
95053acd3b1SBarry Smith 
95153acd3b1SBarry Smith    Input Parameters:
95253acd3b1SBarry Smith +  opt - option name
95353acd3b1SBarry Smith .  text - short string that describes the option
95453acd3b1SBarry Smith .  man - manual page with additional information on option
9550fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
9560fdccdaeSBarry Smith $                 PetscOptionsScalar(..., obj->value,&object->value,...) or
9570fdccdaeSBarry Smith $                 value = defaultvalue
9580fdccdaeSBarry Smith $                 PetscOptionsScalar(..., value,&value,&flg);
9590fdccdaeSBarry Smith $                 if (flg) {
9600fdccdaeSBarry Smith 
96153acd3b1SBarry Smith 
96253acd3b1SBarry Smith    Output Parameter:
96353acd3b1SBarry Smith +  value - the value to return
96453acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
96553acd3b1SBarry Smith 
96695452b02SPatrick Sanan    Notes:
96795452b02SPatrick Sanan     If the user does not supply the option at all value is NOT changed. Thus
9682efd9cb1SBarry Smith           you should ALWAYS initialize value if you access it without first checking if the set flag is true.
9692efd9cb1SBarry Smith 
970989712b9SBarry Smith           The default/currentvalue passed into this routine does not get transferred to the output value variable automatically.
971989712b9SBarry Smith 
97253acd3b1SBarry Smith    Level: beginner
97353acd3b1SBarry Smith 
97495452b02SPatrick Sanan    Notes:
97595452b02SPatrick Sanan     Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
97653acd3b1SBarry Smith 
97789a13869SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
978acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
979acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
98053acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
98153acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
982acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
983a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
98453acd3b1SBarry Smith @*/
9854416b707SBarry Smith PetscErrorCode  PetscOptionsScalar_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscScalar currentvalue,PetscScalar *value,PetscBool  *set)
98653acd3b1SBarry Smith {
98753acd3b1SBarry Smith   PetscErrorCode ierr;
98853acd3b1SBarry Smith 
98953acd3b1SBarry Smith   PetscFunctionBegin;
99053acd3b1SBarry Smith #if !defined(PETSC_USE_COMPLEX)
9910fdccdaeSBarry Smith   ierr = PetscOptionsReal(opt,text,man,currentvalue,value,set);CHKERRQ(ierr);
99253acd3b1SBarry Smith #else
993c5929fdfSBarry Smith   ierr = PetscOptionsGetScalar(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,set);CHKERRQ(ierr);
99453acd3b1SBarry Smith #endif
99553acd3b1SBarry Smith   PetscFunctionReturn(0);
99653acd3b1SBarry Smith }
99753acd3b1SBarry Smith 
99853acd3b1SBarry Smith /*@C
99990d69ab7SBarry 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
100090d69ab7SBarry Smith                       its value is set to false.
100153acd3b1SBarry Smith 
10023f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
100353acd3b1SBarry Smith 
100453acd3b1SBarry Smith    Input Parameters:
100553acd3b1SBarry Smith +  opt - option name
100653acd3b1SBarry Smith .  text - short string that describes the option
100753acd3b1SBarry Smith -  man - manual page with additional information on option
100853acd3b1SBarry Smith 
100953acd3b1SBarry Smith    Output Parameter:
101053acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
101153acd3b1SBarry Smith 
101253acd3b1SBarry Smith    Level: beginner
101353acd3b1SBarry Smith 
101495452b02SPatrick Sanan    Notes:
101595452b02SPatrick Sanan     Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
101653acd3b1SBarry Smith 
101789a13869SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
1018acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
1019acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
102053acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
102153acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1022acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1023a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
102453acd3b1SBarry Smith @*/
10254416b707SBarry Smith PetscErrorCode  PetscOptionsName_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool  *flg)
102653acd3b1SBarry Smith {
102753acd3b1SBarry Smith   PetscErrorCode  ierr;
10284416b707SBarry Smith   PetscOptionItem amsopt;
102953acd3b1SBarry Smith 
103053acd3b1SBarry Smith   PetscFunctionBegin;
1031e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
10324416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1033ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1034a297a907SKarl Rupp 
1035ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
10361ae3d29cSBarry Smith   }
1037c5929fdfSBarry Smith   ierr = PetscOptionsHasName(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,flg);CHKERRQ(ierr);
1038e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1039e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
104053acd3b1SBarry Smith   }
104153acd3b1SBarry Smith   PetscFunctionReturn(0);
104253acd3b1SBarry Smith }
104353acd3b1SBarry Smith 
104453acd3b1SBarry Smith /*@C
1045a264d7a6SBarry Smith      PetscOptionsFList - Puts a list of option values that a single one may be selected from
104653acd3b1SBarry Smith 
10473f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
104853acd3b1SBarry Smith 
104953acd3b1SBarry Smith    Input Parameters:
105053acd3b1SBarry Smith +  opt - option name
105153acd3b1SBarry Smith .  text - short string that describes the option
105253acd3b1SBarry Smith .  man - manual page with additional information on option
105353acd3b1SBarry Smith .  list - the possible choices
10540fdccdaeSBarry Smith .  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with
10550fdccdaeSBarry Smith $                 PetscOptionsFlist(..., obj->value,value,len,&flg);
10560fdccdaeSBarry Smith $                 if (flg) {
10573cc1e11dSBarry Smith -  len - the length of the character array value
105853acd3b1SBarry Smith 
105953acd3b1SBarry Smith    Output Parameter:
106053acd3b1SBarry Smith +  value - the value to return
106153acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
106253acd3b1SBarry Smith 
106353acd3b1SBarry Smith    Level: intermediate
106453acd3b1SBarry Smith 
106595452b02SPatrick Sanan    Notes:
106695452b02SPatrick Sanan     Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
106753acd3b1SBarry Smith 
10682efd9cb1SBarry Smith           If the user does not supply the option at all value is NOT changed. Thus
10692efd9cb1SBarry Smith           you should ALWAYS initialize value if you access it without first checking if the set flag is true.
10702efd9cb1SBarry Smith 
1071989712b9SBarry Smith           The default/currentvalue passed into this routine does not get transferred to the output value variable automatically.
1072989712b9SBarry Smith 
107353acd3b1SBarry Smith    See PetscOptionsEList() for when the choices are given in a string array
107453acd3b1SBarry Smith 
107553acd3b1SBarry Smith    To get a listing of all currently specified options,
107688c29154SBarry Smith     see PetscOptionsView() or PetscOptionsGetAll()
107753acd3b1SBarry Smith 
1078eabe10d7SBarry Smith    Developer Note: This cannot check for invalid selection because of things like MATAIJ that are not included in the list
1079eabe10d7SBarry Smith 
108089a13869SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1081acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
108253acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
108353acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1084acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1085a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsEnum()
108653acd3b1SBarry Smith @*/
10874416b707SBarry 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)
108853acd3b1SBarry Smith {
108953acd3b1SBarry Smith   PetscErrorCode  ierr;
10904416b707SBarry Smith   PetscOptionItem amsopt;
109144ef3d73SBarry Smith   PetscBool       lset;
109253acd3b1SBarry Smith 
109353acd3b1SBarry Smith   PetscFunctionBegin;
10941a1499c8SBarry Smith   if (!PetscOptionsObject->count) {
10954416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,ltext,man,OPTION_FLIST,&amsopt);CHKERRQ(ierr);
109664facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
10970fdccdaeSBarry Smith     ierr = PetscStrdup(currentvalue ? currentvalue : "",(char**)&amsopt->data);CHKERRQ(ierr);
10983cc1e11dSBarry Smith     amsopt->flist = list;
10993cc1e11dSBarry Smith   }
110044ef3d73SBarry Smith   ierr = PetscOptionsGetString(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,len,&lset);CHKERRQ(ierr);
110144ef3d73SBarry Smith   if (set) *set = lset;
11021a1499c8SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
110344ef3d73SBarry Smith     ierr = PetscFunctionListPrintTypes(PetscOptionsObject->comm,stdout,PetscOptionsObject->prefix,opt,ltext,man,list,currentvalue,lset && value ? value : currentvalue);CHKERRQ(ierr);CHKERRQ(ierr);
110453acd3b1SBarry Smith   }
110553acd3b1SBarry Smith   PetscFunctionReturn(0);
110653acd3b1SBarry Smith }
110753acd3b1SBarry Smith 
110853acd3b1SBarry Smith /*@C
110953acd3b1SBarry Smith      PetscOptionsEList - Puts a list of option values that a single one may be selected from
111053acd3b1SBarry Smith 
11113f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
111253acd3b1SBarry Smith 
111353acd3b1SBarry Smith    Input Parameters:
111453acd3b1SBarry Smith +  opt - option name
111553acd3b1SBarry Smith .  ltext - short string that describes the option
111653acd3b1SBarry Smith .  man - manual page with additional information on option
1117a264d7a6SBarry Smith .  list - the possible choices (one of these must be selected, anything else is invalid)
111853acd3b1SBarry Smith .  ntext - number of choices
11190fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with
11200fdccdaeSBarry Smith $                 PetscOptionsElist(..., obj->value,&value,&flg);
11210fdccdaeSBarry Smith $                 if (flg) {
11220fdccdaeSBarry Smith 
112353acd3b1SBarry Smith 
112453acd3b1SBarry Smith    Output Parameter:
112553acd3b1SBarry Smith +  value - the index of the value to return
112653acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
112753acd3b1SBarry Smith 
112853acd3b1SBarry Smith    Level: intermediate
112953acd3b1SBarry Smith 
113095452b02SPatrick Sanan    Notes:
113195452b02SPatrick Sanan     Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
113253acd3b1SBarry Smith 
11332efd9cb1SBarry Smith          If the user does not supply the option at all value is NOT changed. Thus
11342efd9cb1SBarry Smith           you should ALWAYS initialize value if you access it without first checking if the set flag is true.
11352efd9cb1SBarry Smith 
1136a264d7a6SBarry Smith    See PetscOptionsFList() for when the choices are given in a PetscFunctionList()
113753acd3b1SBarry Smith 
113889a13869SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1139acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
114053acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
114153acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1142acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1143a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEnum()
114453acd3b1SBarry Smith @*/
11454416b707SBarry 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)
114653acd3b1SBarry Smith {
114753acd3b1SBarry Smith   PetscErrorCode  ierr;
114853acd3b1SBarry Smith   PetscInt        i;
11494416b707SBarry Smith   PetscOptionItem amsopt;
115044ef3d73SBarry Smith   PetscBool       lset;
115153acd3b1SBarry Smith 
115253acd3b1SBarry Smith   PetscFunctionBegin;
11531a1499c8SBarry Smith   if (!PetscOptionsObject->count) {
11544416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,ltext,man,OPTION_ELIST,&amsopt);CHKERRQ(ierr);
115564facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
11560fdccdaeSBarry Smith     ierr = PetscStrdup(currentvalue ? currentvalue : "",(char**)&amsopt->data);CHKERRQ(ierr);
11576991f827SBarry Smith     ierr = PetscStrNArrayallocpy(ntext,list,(char***)&amsopt->list);CHKERRQ(ierr);
11581ae3d29cSBarry Smith     amsopt->nlist = ntext;
11591ae3d29cSBarry Smith   }
116044ef3d73SBarry Smith   ierr = PetscOptionsGetEList(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,list,ntext,value,&lset);CHKERRQ(ierr);
116144ef3d73SBarry Smith   if (set) *set = lset;
11621a1499c8SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1163*0e4c290aSBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <now %s : formerly %s> %s (choose one of)",PetscOptionsObject->prefix?PetscOptionsObject->prefix:"",opt+1,lset && value ? list[*value] : currentvalue,currentvalue,ltext);CHKERRQ(ierr);
116453acd3b1SBarry Smith     for (i=0; i<ntext; i++) {
1165e55864a3SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," %s",list[i]);CHKERRQ(ierr);
116653acd3b1SBarry Smith     }
1167e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," (%s)\n",ManSection(man));CHKERRQ(ierr);
116853acd3b1SBarry Smith   }
116953acd3b1SBarry Smith   PetscFunctionReturn(0);
117053acd3b1SBarry Smith }
117153acd3b1SBarry Smith 
117253acd3b1SBarry Smith /*@C
1173acfcf0e5SJed Brown      PetscOptionsBoolGroupBegin - First in a series of logical queries on the options database for
1174d5649816SBarry Smith        which at most a single value can be true.
117553acd3b1SBarry Smith 
11763f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
117753acd3b1SBarry Smith 
117853acd3b1SBarry Smith    Input Parameters:
117953acd3b1SBarry Smith +  opt - option name
118053acd3b1SBarry Smith .  text - short string that describes the option
118153acd3b1SBarry Smith -  man - manual page with additional information on option
118253acd3b1SBarry Smith 
118353acd3b1SBarry Smith    Output Parameter:
118453acd3b1SBarry Smith .  flg - whether that option was set or not
118553acd3b1SBarry Smith 
118653acd3b1SBarry Smith    Level: intermediate
118753acd3b1SBarry Smith 
118895452b02SPatrick Sanan    Notes:
118995452b02SPatrick Sanan     Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
119053acd3b1SBarry Smith 
1191acfcf0e5SJed Brown    Must be followed by 0 or more PetscOptionsBoolGroup()s and PetscOptionsBoolGroupEnd()
119253acd3b1SBarry Smith 
119389a13869SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1194acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
119553acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
119653acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1197acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1198a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
119953acd3b1SBarry Smith @*/
12004416b707SBarry Smith PetscErrorCode  PetscOptionsBoolGroupBegin_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool  *flg)
120153acd3b1SBarry Smith {
120253acd3b1SBarry Smith   PetscErrorCode  ierr;
12034416b707SBarry Smith   PetscOptionItem amsopt;
120453acd3b1SBarry Smith 
120553acd3b1SBarry Smith   PetscFunctionBegin;
1206e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
12074416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1208ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1209a297a907SKarl Rupp 
1210ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
12111ae3d29cSBarry Smith   }
121268b16fdaSBarry Smith   *flg = PETSC_FALSE;
1213c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,flg,NULL);CHKERRQ(ierr);
1214e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1215e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  Pick at most one of -------------\n");CHKERRQ(ierr);
1216e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"    -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
121753acd3b1SBarry Smith   }
121853acd3b1SBarry Smith   PetscFunctionReturn(0);
121953acd3b1SBarry Smith }
122053acd3b1SBarry Smith 
122153acd3b1SBarry Smith /*@C
1222acfcf0e5SJed Brown      PetscOptionsBoolGroup - One in a series of logical queries on the options database for
1223d5649816SBarry Smith        which at most a single value can be true.
122453acd3b1SBarry Smith 
12253f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
122653acd3b1SBarry Smith 
122753acd3b1SBarry Smith    Input Parameters:
122853acd3b1SBarry Smith +  opt - option name
122953acd3b1SBarry Smith .  text - short string that describes the option
123053acd3b1SBarry Smith -  man - manual page with additional information on option
123153acd3b1SBarry Smith 
123253acd3b1SBarry Smith    Output Parameter:
123353acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
123453acd3b1SBarry Smith 
123553acd3b1SBarry Smith    Level: intermediate
123653acd3b1SBarry Smith 
123795452b02SPatrick Sanan    Notes:
123895452b02SPatrick Sanan     Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
123953acd3b1SBarry Smith 
1240acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin() and preceded a PetscOptionsBoolGroupEnd()
124153acd3b1SBarry Smith 
124289a13869SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1243acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
124453acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
124553acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1246acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1247a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
124853acd3b1SBarry Smith @*/
12494416b707SBarry Smith PetscErrorCode  PetscOptionsBoolGroup_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool  *flg)
125053acd3b1SBarry Smith {
125153acd3b1SBarry Smith   PetscErrorCode  ierr;
12524416b707SBarry Smith   PetscOptionItem amsopt;
125353acd3b1SBarry Smith 
125453acd3b1SBarry Smith   PetscFunctionBegin;
1255e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
12564416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1257ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1258a297a907SKarl Rupp 
1259ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
12601ae3d29cSBarry Smith   }
126117326d04SJed Brown   *flg = PETSC_FALSE;
1262c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,flg,NULL);CHKERRQ(ierr);
1263e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1264e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"    -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
126553acd3b1SBarry Smith   }
126653acd3b1SBarry Smith   PetscFunctionReturn(0);
126753acd3b1SBarry Smith }
126853acd3b1SBarry Smith 
126953acd3b1SBarry Smith /*@C
1270acfcf0e5SJed Brown      PetscOptionsBoolGroupEnd - Last in a series of logical queries on the options database for
1271d5649816SBarry Smith        which at most a single value can be true.
127253acd3b1SBarry Smith 
12733f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
127453acd3b1SBarry Smith 
127553acd3b1SBarry Smith    Input Parameters:
127653acd3b1SBarry Smith +  opt - option name
127753acd3b1SBarry Smith .  text - short string that describes the option
127853acd3b1SBarry Smith -  man - manual page with additional information on option
127953acd3b1SBarry Smith 
128053acd3b1SBarry Smith    Output Parameter:
128153acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
128253acd3b1SBarry Smith 
128353acd3b1SBarry Smith    Level: intermediate
128453acd3b1SBarry Smith 
128595452b02SPatrick Sanan    Notes:
128695452b02SPatrick Sanan     Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
128753acd3b1SBarry Smith 
1288acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin()
128953acd3b1SBarry Smith 
129089a13869SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1291acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
129253acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
129353acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1294acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1295a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
129653acd3b1SBarry Smith @*/
12974416b707SBarry Smith PetscErrorCode  PetscOptionsBoolGroupEnd_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool  *flg)
129853acd3b1SBarry Smith {
129953acd3b1SBarry Smith   PetscErrorCode  ierr;
13004416b707SBarry Smith   PetscOptionItem amsopt;
130153acd3b1SBarry Smith 
130253acd3b1SBarry Smith   PetscFunctionBegin;
1303e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
13044416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1305ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1306a297a907SKarl Rupp 
1307ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
13081ae3d29cSBarry Smith   }
130917326d04SJed Brown   *flg = PETSC_FALSE;
1310c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,flg,NULL);CHKERRQ(ierr);
1311e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1312e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"    -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
131353acd3b1SBarry Smith   }
131453acd3b1SBarry Smith   PetscFunctionReturn(0);
131553acd3b1SBarry Smith }
131653acd3b1SBarry Smith 
131753acd3b1SBarry Smith /*@C
1318acfcf0e5SJed Brown    PetscOptionsBool - Determines if a particular option is in the database with a true or false
131953acd3b1SBarry Smith 
13203f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
132153acd3b1SBarry Smith 
132253acd3b1SBarry Smith    Input Parameters:
132353acd3b1SBarry Smith +  opt - option name
132453acd3b1SBarry Smith .  text - short string that describes the option
1325868c398cSBarry Smith .  man - manual page with additional information on option
132694ae4db5SBarry Smith -  currentvalue - the current value
132753acd3b1SBarry Smith 
132853acd3b1SBarry Smith    Output Parameter:
132953acd3b1SBarry Smith .  flg - PETSC_TRUE or PETSC_FALSE
133053acd3b1SBarry Smith .  set - PETSC_TRUE if found, else PETSC_FALSE
133153acd3b1SBarry Smith 
13322efd9cb1SBarry Smith    Notes:
13332efd9cb1SBarry Smith        TRUE, true, YES, yes, nostring, and 1 all translate to PETSC_TRUE
13342efd9cb1SBarry Smith        FALSE, false, NO, no, and 0 all translate to PETSC_FALSE
13352efd9cb1SBarry Smith 
13362efd9cb1SBarry 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
13372efd9cb1SBarry Smith      is equivalent to -requested_bool true
13382efd9cb1SBarry Smith 
13392efd9cb1SBarry Smith        If the user does not supply the option at all flg is NOT changed. Thus
13402efd9cb1SBarry Smith      you should ALWAYS initialize the flg if you access it without first checking if the set flag is true.
13412efd9cb1SBarry Smith 
134253acd3b1SBarry Smith    Level: beginner
134353acd3b1SBarry Smith 
134495452b02SPatrick Sanan    Notes:
134595452b02SPatrick Sanan     Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
134653acd3b1SBarry Smith 
134789a13869SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
1348acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
1349acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
135053acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
135153acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1352acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1353a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
135453acd3b1SBarry Smith @*/
13554416b707SBarry Smith PetscErrorCode  PetscOptionsBool_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool currentvalue,PetscBool  *flg,PetscBool  *set)
135653acd3b1SBarry Smith {
135753acd3b1SBarry Smith   PetscErrorCode  ierr;
1358ace3abfcSBarry Smith   PetscBool       iset;
13594416b707SBarry Smith   PetscOptionItem amsopt;
136053acd3b1SBarry Smith 
136153acd3b1SBarry Smith   PetscFunctionBegin;
1362e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
13634416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1364ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1365a297a907SKarl Rupp 
136694ae4db5SBarry Smith     *(PetscBool*)amsopt->data = currentvalue;
1367af6d86caSBarry Smith   }
1368c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,flg,&iset);CHKERRQ(ierr);
136953acd3b1SBarry Smith   if (set) *set = iset;
13701a1499c8SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
137144ef3d73SBarry Smith     const char *v = PetscBools[currentvalue], *vn = PetscBools[iset && flg ? *flg : currentvalue];
137244ef3d73SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s: <%s : %s> %s (%s)\n",PetscOptionsObject->prefix?PetscOptionsObject->prefix:"",opt+1,v,vn,text,ManSection(man));CHKERRQ(ierr);
137353acd3b1SBarry Smith   }
137453acd3b1SBarry Smith   PetscFunctionReturn(0);
137553acd3b1SBarry Smith }
137653acd3b1SBarry Smith 
137753acd3b1SBarry Smith /*@C
137853acd3b1SBarry Smith    PetscOptionsRealArray - Gets an array of double values for a particular
137953acd3b1SBarry Smith    option in the database. The values must be separated with commas with
138053acd3b1SBarry Smith    no intervening spaces.
138153acd3b1SBarry Smith 
13823f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
138353acd3b1SBarry Smith 
138453acd3b1SBarry Smith    Input Parameters:
138553acd3b1SBarry Smith +  opt - the option one is seeking
138653acd3b1SBarry Smith .  text - short string describing option
138753acd3b1SBarry Smith .  man - manual page for option
138853acd3b1SBarry Smith -  nmax - maximum number of values
138953acd3b1SBarry Smith 
139053acd3b1SBarry Smith    Output Parameter:
139153acd3b1SBarry Smith +  value - location to copy values
139253acd3b1SBarry Smith .  nmax - actual number of values found
139353acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
139453acd3b1SBarry Smith 
139553acd3b1SBarry Smith    Level: beginner
139653acd3b1SBarry Smith 
139753acd3b1SBarry Smith    Notes:
139853acd3b1SBarry Smith    The user should pass in an array of doubles
139953acd3b1SBarry Smith 
140053acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
140153acd3b1SBarry Smith 
140289a13869SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1403acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
140453acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
140553acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1406acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1407a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
140853acd3b1SBarry Smith @*/
14094416b707SBarry Smith PetscErrorCode PetscOptionsRealArray_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscReal value[],PetscInt *n,PetscBool  *set)
141053acd3b1SBarry Smith {
141153acd3b1SBarry Smith   PetscErrorCode  ierr;
141253acd3b1SBarry Smith   PetscInt        i;
14134416b707SBarry Smith   PetscOptionItem amsopt;
141453acd3b1SBarry Smith 
141553acd3b1SBarry Smith   PetscFunctionBegin;
1416e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
1417e26ddf31SBarry Smith     PetscReal *vals;
1418e26ddf31SBarry Smith 
14194416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_REAL_ARRAY,&amsopt);CHKERRQ(ierr);
1420e55864a3SBarry Smith     ierr = PetscMalloc((*n)*sizeof(PetscReal),&amsopt->data);CHKERRQ(ierr);
1421e26ddf31SBarry Smith     vals = (PetscReal*)amsopt->data;
1422e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1423e26ddf31SBarry Smith     amsopt->arraylength = *n;
1424e26ddf31SBarry Smith   }
1425c5929fdfSBarry Smith   ierr = PetscOptionsGetRealArray(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,n,set);CHKERRQ(ierr);
1426e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1427a519f713SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%g",PetscOptionsObject->prefix?PetscOptionsObject->prefix:"",opt+1,(double)value[0]);CHKERRQ(ierr);
142853acd3b1SBarry Smith     for (i=1; i<*n; i++) {
1429e55864a3SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%g",(double)value[i]);CHKERRQ(ierr);
143053acd3b1SBarry Smith     }
1431e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
143253acd3b1SBarry Smith   }
143353acd3b1SBarry Smith   PetscFunctionReturn(0);
143453acd3b1SBarry Smith }
143553acd3b1SBarry Smith 
1436050cccc3SHong Zhang /*@C
1437050cccc3SHong Zhang    PetscOptionsScalarArray - Gets an array of Scalar values for a particular
1438050cccc3SHong Zhang    option in the database. The values must be separated with commas with
1439050cccc3SHong Zhang    no intervening spaces.
1440050cccc3SHong Zhang 
1441050cccc3SHong Zhang    Logically Collective on the communicator passed in PetscOptionsBegin()
1442050cccc3SHong Zhang 
1443050cccc3SHong Zhang    Input Parameters:
1444050cccc3SHong Zhang +  opt - the option one is seeking
1445050cccc3SHong Zhang .  text - short string describing option
1446050cccc3SHong Zhang .  man - manual page for option
1447050cccc3SHong Zhang -  nmax - maximum number of values
1448050cccc3SHong Zhang 
1449050cccc3SHong Zhang    Output Parameter:
1450050cccc3SHong Zhang +  value - location to copy values
1451050cccc3SHong Zhang .  nmax - actual number of values found
1452050cccc3SHong Zhang -  set - PETSC_TRUE if found, else PETSC_FALSE
1453050cccc3SHong Zhang 
1454050cccc3SHong Zhang    Level: beginner
1455050cccc3SHong Zhang 
1456050cccc3SHong Zhang    Notes:
1457050cccc3SHong Zhang    The user should pass in an array of doubles
1458050cccc3SHong Zhang 
1459050cccc3SHong Zhang    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
1460050cccc3SHong Zhang 
146189a13869SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1462050cccc3SHong Zhang           PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
1463050cccc3SHong Zhang           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1464050cccc3SHong Zhang           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1465050cccc3SHong Zhang           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1466050cccc3SHong Zhang           PetscOptionsFList(), PetscOptionsEList()
1467050cccc3SHong Zhang @*/
14684416b707SBarry Smith PetscErrorCode PetscOptionsScalarArray_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscScalar value[],PetscInt *n,PetscBool  *set)
1469050cccc3SHong Zhang {
1470050cccc3SHong Zhang   PetscErrorCode  ierr;
1471050cccc3SHong Zhang   PetscInt        i;
14724416b707SBarry Smith   PetscOptionItem amsopt;
1473050cccc3SHong Zhang 
1474050cccc3SHong Zhang   PetscFunctionBegin;
1475050cccc3SHong Zhang   if (!PetscOptionsObject->count) {
1476050cccc3SHong Zhang     PetscScalar *vals;
1477050cccc3SHong Zhang 
14784416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_SCALAR_ARRAY,&amsopt);CHKERRQ(ierr);
1479050cccc3SHong Zhang     ierr = PetscMalloc((*n)*sizeof(PetscScalar),&amsopt->data);CHKERRQ(ierr);
1480050cccc3SHong Zhang     vals = (PetscScalar*)amsopt->data;
1481050cccc3SHong Zhang     for (i=0; i<*n; i++) vals[i] = value[i];
1482050cccc3SHong Zhang     amsopt->arraylength = *n;
1483050cccc3SHong Zhang   }
1484c5929fdfSBarry Smith   ierr = PetscOptionsGetScalarArray(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,n,set);CHKERRQ(ierr);
1485050cccc3SHong Zhang   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
148695f3a755SJose 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);
1487050cccc3SHong Zhang     for (i=1; i<*n; i++) {
148895f3a755SJose E. Roman       ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%g+%gi",(double)PetscRealPart(value[i]),(double)PetscImaginaryPart(value[i]));CHKERRQ(ierr);
1489050cccc3SHong Zhang     }
1490050cccc3SHong Zhang     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
1491050cccc3SHong Zhang   }
1492050cccc3SHong Zhang   PetscFunctionReturn(0);
1493050cccc3SHong Zhang }
149453acd3b1SBarry Smith 
149553acd3b1SBarry Smith /*@C
149653acd3b1SBarry Smith    PetscOptionsIntArray - Gets an array of integers for a particular
1497b32a342fSShri Abhyankar    option in the database.
149853acd3b1SBarry Smith 
14993f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
150053acd3b1SBarry Smith 
150153acd3b1SBarry Smith    Input Parameters:
150253acd3b1SBarry Smith +  opt - the option one is seeking
150353acd3b1SBarry Smith .  text - short string describing option
150453acd3b1SBarry Smith .  man - manual page for option
1505f8a50e2bSBarry Smith -  n - maximum number of values
150653acd3b1SBarry Smith 
150753acd3b1SBarry Smith    Output Parameter:
150853acd3b1SBarry Smith +  value - location to copy values
1509f8a50e2bSBarry Smith .  n - actual number of values found
151053acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
151153acd3b1SBarry Smith 
151253acd3b1SBarry Smith    Level: beginner
151353acd3b1SBarry Smith 
151453acd3b1SBarry Smith    Notes:
1515b32a342fSShri Abhyankar    The array can be passed as
1516bebe2cf6SSatish Balay    a comma separated list:                                 0,1,2,3,4,5,6,7
15170fd488f5SShri Abhyankar    a range (start-end+1):                                  0-8
15180fd488f5SShri Abhyankar    a range with given increment (start-end+1:inc):         0-7:2
1519bebe2cf6SSatish Balay    a combination of values and ranges separated by commas: 0,1-8,8-15:2
1520b32a342fSShri Abhyankar 
1521b32a342fSShri Abhyankar    There must be no intervening spaces between the values.
152253acd3b1SBarry Smith 
152353acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
152453acd3b1SBarry Smith 
152589a13869SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1526acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
152753acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
152853acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1529acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1530a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsRealArray()
153153acd3b1SBarry Smith @*/
15324416b707SBarry Smith PetscErrorCode  PetscOptionsIntArray_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscInt value[],PetscInt *n,PetscBool  *set)
153353acd3b1SBarry Smith {
153453acd3b1SBarry Smith   PetscErrorCode ierr;
153553acd3b1SBarry Smith   PetscInt        i;
15364416b707SBarry Smith   PetscOptionItem amsopt;
153753acd3b1SBarry Smith 
153853acd3b1SBarry Smith   PetscFunctionBegin;
1539e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
1540e26ddf31SBarry Smith     PetscInt *vals;
1541e26ddf31SBarry Smith 
15424416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_INT_ARRAY,&amsopt);CHKERRQ(ierr);
1543854ce69bSBarry Smith     ierr = PetscMalloc1(*n,(PetscInt**)&amsopt->data);CHKERRQ(ierr);
1544e26ddf31SBarry Smith     vals = (PetscInt*)amsopt->data;
1545e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1546e26ddf31SBarry Smith     amsopt->arraylength = *n;
1547e26ddf31SBarry Smith   }
1548c5929fdfSBarry Smith   ierr = PetscOptionsGetIntArray(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,n,set);CHKERRQ(ierr);
1549e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1550e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%d",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,value[0]);CHKERRQ(ierr);
155153acd3b1SBarry Smith     for (i=1; i<*n; i++) {
1552e55864a3SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%d",value[i]);CHKERRQ(ierr);
155353acd3b1SBarry Smith     }
1554e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
155553acd3b1SBarry Smith   }
155653acd3b1SBarry Smith   PetscFunctionReturn(0);
155753acd3b1SBarry Smith }
155853acd3b1SBarry Smith 
155953acd3b1SBarry Smith /*@C
156053acd3b1SBarry Smith    PetscOptionsStringArray - Gets an array of string values for a particular
156153acd3b1SBarry Smith    option in the database. The values must be separated with commas with
156253acd3b1SBarry Smith    no intervening spaces.
156353acd3b1SBarry Smith 
15643f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
156553acd3b1SBarry Smith 
156653acd3b1SBarry Smith    Input Parameters:
156753acd3b1SBarry Smith +  opt - the option one is seeking
156853acd3b1SBarry Smith .  text - short string describing option
156953acd3b1SBarry Smith .  man - manual page for option
157053acd3b1SBarry Smith -  nmax - maximum number of strings
157153acd3b1SBarry Smith 
157253acd3b1SBarry Smith    Output Parameter:
157353acd3b1SBarry Smith +  value - location to copy strings
157453acd3b1SBarry Smith .  nmax - actual number of strings found
157553acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
157653acd3b1SBarry Smith 
157753acd3b1SBarry Smith    Level: beginner
157853acd3b1SBarry Smith 
157953acd3b1SBarry Smith    Notes:
158053acd3b1SBarry Smith    The user should pass in an array of pointers to char, to hold all the
158153acd3b1SBarry Smith    strings returned by this function.
158253acd3b1SBarry Smith 
158353acd3b1SBarry Smith    The user is responsible for deallocating the strings that are
158453acd3b1SBarry Smith    returned. The Fortran interface for this routine is not supported.
158553acd3b1SBarry Smith 
158653acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
158753acd3b1SBarry Smith 
158889a13869SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1589acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
159053acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
159153acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1592acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1593a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
159453acd3b1SBarry Smith @*/
15954416b707SBarry Smith PetscErrorCode  PetscOptionsStringArray_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],char *value[],PetscInt *nmax,PetscBool  *set)
159653acd3b1SBarry Smith {
159753acd3b1SBarry Smith   PetscErrorCode  ierr;
15984416b707SBarry Smith   PetscOptionItem amsopt;
159953acd3b1SBarry Smith 
160053acd3b1SBarry Smith   PetscFunctionBegin;
1601e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
16024416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_STRING_ARRAY,&amsopt);CHKERRQ(ierr);
1603854ce69bSBarry Smith     ierr = PetscMalloc1(*nmax,(char**)&amsopt->data);CHKERRQ(ierr);
1604a297a907SKarl Rupp 
16051ae3d29cSBarry Smith     amsopt->arraylength = *nmax;
16061ae3d29cSBarry Smith   }
1607c5929fdfSBarry Smith   ierr = PetscOptionsGetStringArray(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,nmax,set);CHKERRQ(ierr);
1608e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1609e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <string1,string2,...>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
161053acd3b1SBarry Smith   }
161153acd3b1SBarry Smith   PetscFunctionReturn(0);
161253acd3b1SBarry Smith }
161353acd3b1SBarry Smith 
1614e2446a98SMatthew Knepley /*@C
1615acfcf0e5SJed Brown    PetscOptionsBoolArray - Gets an array of logical values (true or false) for a particular
1616e2446a98SMatthew Knepley    option in the database. The values must be separated with commas with
1617e2446a98SMatthew Knepley    no intervening spaces.
1618e2446a98SMatthew Knepley 
16193f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
1620e2446a98SMatthew Knepley 
1621e2446a98SMatthew Knepley    Input Parameters:
1622e2446a98SMatthew Knepley +  opt - the option one is seeking
1623e2446a98SMatthew Knepley .  text - short string describing option
1624e2446a98SMatthew Knepley .  man - manual page for option
1625e2446a98SMatthew Knepley -  nmax - maximum number of values
1626e2446a98SMatthew Knepley 
1627e2446a98SMatthew Knepley    Output Parameter:
1628e2446a98SMatthew Knepley +  value - location to copy values
1629e2446a98SMatthew Knepley .  nmax - actual number of values found
1630e2446a98SMatthew Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
1631e2446a98SMatthew Knepley 
1632e2446a98SMatthew Knepley    Level: beginner
1633e2446a98SMatthew Knepley 
1634e2446a98SMatthew Knepley    Notes:
1635e2446a98SMatthew Knepley    The user should pass in an array of doubles
1636e2446a98SMatthew Knepley 
1637e2446a98SMatthew Knepley    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
1638e2446a98SMatthew Knepley 
163989a13869SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1640acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
1641e2446a98SMatthew Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1642e2446a98SMatthew Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1643acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1644a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
1645e2446a98SMatthew Knepley @*/
16464416b707SBarry Smith PetscErrorCode  PetscOptionsBoolArray_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool value[],PetscInt *n,PetscBool *set)
1647e2446a98SMatthew Knepley {
1648e2446a98SMatthew Knepley   PetscErrorCode   ierr;
1649e2446a98SMatthew Knepley   PetscInt         i;
16504416b707SBarry Smith   PetscOptionItem  amsopt;
1651e2446a98SMatthew Knepley 
1652e2446a98SMatthew Knepley   PetscFunctionBegin;
1653e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
1654ace3abfcSBarry Smith     PetscBool *vals;
16551ae3d29cSBarry Smith 
16564416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL_ARRAY,&amsopt);CHKERRQ(ierr);
16571a1499c8SBarry Smith     ierr = PetscMalloc1(*n,(PetscBool**)&amsopt->data);CHKERRQ(ierr);
1658ace3abfcSBarry Smith     vals = (PetscBool*)amsopt->data;
16591ae3d29cSBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
16601ae3d29cSBarry Smith     amsopt->arraylength = *n;
16611ae3d29cSBarry Smith   }
1662c5929fdfSBarry Smith   ierr = PetscOptionsGetBoolArray(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,n,set);CHKERRQ(ierr);
1663e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1664e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%d",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,value[0]);CHKERRQ(ierr);
1665e2446a98SMatthew Knepley     for (i=1; i<*n; i++) {
1666e55864a3SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%d",value[i]);CHKERRQ(ierr);
1667e2446a98SMatthew Knepley     }
1668e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
1669e2446a98SMatthew Knepley   }
1670e2446a98SMatthew Knepley   PetscFunctionReturn(0);
1671e2446a98SMatthew Knepley }
1672e2446a98SMatthew Knepley 
16738cc676e6SMatthew G Knepley /*@C
1674d1da0b69SBarry Smith    PetscOptionsViewer - Gets a viewer appropriate for the type indicated by the user
16758cc676e6SMatthew G Knepley 
16768cc676e6SMatthew G Knepley    Logically Collective on the communicator passed in PetscOptionsBegin()
16778cc676e6SMatthew G Knepley 
16788cc676e6SMatthew G Knepley    Input Parameters:
16798cc676e6SMatthew G Knepley +  opt - option name
16808cc676e6SMatthew G Knepley .  text - short string that describes the option
16818cc676e6SMatthew G Knepley -  man - manual page with additional information on option
16828cc676e6SMatthew G Knepley 
16838cc676e6SMatthew G Knepley    Output Parameter:
16848cc676e6SMatthew G Knepley +  viewer - the viewer
16858cc676e6SMatthew G Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
16868cc676e6SMatthew G Knepley 
16878cc676e6SMatthew G Knepley    Level: beginner
16888cc676e6SMatthew G Knepley 
168995452b02SPatrick Sanan    Notes:
169095452b02SPatrick Sanan     Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
16918cc676e6SMatthew G Knepley 
16925a7113b9SPatrick Sanan    See PetscOptionsGetViewer() for the format of the supplied viewer and its options
16938cc676e6SMatthew G Knepley 
169489a13869SBarry Smith .seealso: PetscOptionsGetViewer(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
16958cc676e6SMatthew G Knepley           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
16968cc676e6SMatthew G Knepley           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
16978cc676e6SMatthew G Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
16988cc676e6SMatthew G Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
16998cc676e6SMatthew G Knepley           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1700a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
17018cc676e6SMatthew G Knepley @*/
17024416b707SBarry Smith PetscErrorCode  PetscOptionsViewer_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscViewer *viewer,PetscViewerFormat *format,PetscBool  *set)
17038cc676e6SMatthew G Knepley {
17048cc676e6SMatthew G Knepley   PetscErrorCode  ierr;
17054416b707SBarry Smith   PetscOptionItem amsopt;
17068cc676e6SMatthew G Knepley 
17078cc676e6SMatthew G Knepley   PetscFunctionBegin;
17081a1499c8SBarry Smith   if (!PetscOptionsObject->count) {
17094416b707SBarry Smith     ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
171064facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
17115b02f95dSBarry Smith     ierr = PetscStrdup("",(char**)&amsopt->data);CHKERRQ(ierr);
17128cc676e6SMatthew G Knepley   }
171316413a6aSBarry Smith   ierr = PetscOptionsGetViewer(PetscOptionsObject->comm,PetscOptionsObject->options,PetscOptionsObject->prefix,opt,viewer,format,set);CHKERRQ(ierr);
1714e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1715e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%s>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,"",text,ManSection(man));CHKERRQ(ierr);
17168cc676e6SMatthew G Knepley   }
17178cc676e6SMatthew G Knepley   PetscFunctionReturn(0);
17188cc676e6SMatthew G Knepley }
17198cc676e6SMatthew G Knepley 
172053acd3b1SBarry Smith 
172153acd3b1SBarry Smith /*@C
1722b52f573bSBarry Smith      PetscOptionsHead - Puts a heading before listing any more published options. Used, for example,
172353acd3b1SBarry Smith             in KSPSetFromOptions_GMRES().
172453acd3b1SBarry Smith 
17253f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
172653acd3b1SBarry Smith 
172753acd3b1SBarry Smith    Input Parameter:
172853acd3b1SBarry Smith .   head - the heading text
172953acd3b1SBarry Smith 
173053acd3b1SBarry Smith 
173153acd3b1SBarry Smith    Level: intermediate
173253acd3b1SBarry Smith 
173395452b02SPatrick Sanan    Notes:
173495452b02SPatrick Sanan     Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
173553acd3b1SBarry Smith 
1736b52f573bSBarry Smith           Can be followed by a call to PetscOptionsTail() in the same function.
173753acd3b1SBarry Smith 
173889a13869SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1739acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
174053acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
174153acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1742acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1743a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
174453acd3b1SBarry Smith @*/
17454416b707SBarry Smith PetscErrorCode  PetscOptionsHead(PetscOptionItems *PetscOptionsObject,const char head[])
174653acd3b1SBarry Smith {
174753acd3b1SBarry Smith   PetscErrorCode ierr;
174853acd3b1SBarry Smith 
174953acd3b1SBarry Smith   PetscFunctionBegin;
1750e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1751e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  %s\n",head);CHKERRQ(ierr);
175253acd3b1SBarry Smith   }
175353acd3b1SBarry Smith   PetscFunctionReturn(0);
175453acd3b1SBarry Smith }
175553acd3b1SBarry Smith 
175653acd3b1SBarry Smith 
175753acd3b1SBarry Smith 
175853acd3b1SBarry Smith 
175953acd3b1SBarry Smith 
176053acd3b1SBarry Smith 
1761