3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
9 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
10 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
11 * Copyright (c) 2001-2009, The GROMACS development team,
12 * check out http://www.gromacs.org for more information.
14 * This program is free software; you can redistribute it and/or
15 * modify it under the terms of the GNU General Public License
16 * as published by the Free Software Foundation; either version 2
17 * of the License, or (at your option) any later version.
19 * If you want to redistribute modifications, please consider that
20 * scientific software is very special. Version control is crucial -
21 * bugs must be traceable. We will be happy to consider code for
22 * inclusion in the official distribution, but derived work must not
23 * be called official GROMACS. Details are found in the README & COPYING
24 * files - if they are missing, get the official version at www.gromacs.org.
26 * To help us fund GROMACS development, we humbly ask that you cite
27 * the papers on the package - you can find them in the top README file.
29 * For more info, check our website at http://www.gromacs.org
33 * Implementation of functions in selparam.h.
44 #include <selmethod.h>
47 #include "parsetree.h"
52 * \param[in] name Name of the parameter to search.
53 * \param[in] nparam Number of parameters in the \p param array.
54 * \param[in] param Parameter array to search.
55 * \returns Pointer to the parameter in the \p param
56 * or NULL if no parameter with name \p name was found.
58 * The comparison is case-sensitive.
61 gmx_ana_selparam_find(const char *name
, int nparam
, gmx_ana_selparam_t
*param
)
69 /* Find the first non-null parameter */
71 while (i
< nparam
&& param
[i
].name
== NULL
)
75 /* Process the special case of a NULL parameter */
78 return (i
== 0) ? NULL
: ¶m
[i
-1];
80 for ( ; i
< nparam
; ++i
)
82 if (!strcmp(param
[i
].name
, name
))
86 /* Check for 'no' prefix on gmx_boolean parameters */
87 if (param
[i
].val
.type
== NO_VALUE
88 && strlen(name
) > 2 && name
[0] == 'n' && name
[1] == 'o'
89 && !strcmp(param
[i
].name
, name
+2))
98 * Does a type conversion on a \c t_selexpr_value.
100 * \param[in,out] value Value to convert.
101 * \param[in] type Type to convert to.
102 * \param[in] scanner Scanner data structure.
103 * \returns 0 on success, a non-zero value on error.
106 convert_value(t_selexpr_value
*value
, e_selvalue_t type
, void *scanner
)
108 if (value
->type
== type
|| type
== NO_VALUE
)
114 /* Conversion from atom selection to position using default
115 * reference positions. */
116 if (value
->type
== GROUP_VALUE
&& type
== POS_VALUE
)
119 _gmx_sel_init_position(value
->u
.expr
, NULL
, scanner
);
120 if (value
->u
.expr
== NULL
)
131 /* Integers to floating point are easy */
132 if (value
->type
== INT_VALUE
&& type
== REAL_VALUE
)
134 value
->u
.r
.r1
= (real
)value
->u
.i
.i1
;
135 value
->u
.r
.r2
= (real
)value
->u
.i
.i2
;
139 /* Reals that are integer-valued can also be converted */
140 if (value
->type
== REAL_VALUE
&& type
== INT_VALUE
141 && gmx_within_tol(value
->u
.r
.r1
, (int)value
->u
.r
.r1
, GMX_REAL_EPS
)
142 && gmx_within_tol(value
->u
.r
.r2
, (int)value
->u
.r
.r2
, GMX_REAL_EPS
))
144 value
->u
.i
.i1
= (int)value
->u
.r
.r1
;
145 value
->u
.i
.i2
= (int)value
->u
.r
.r2
;
154 * Does a type conversion on a list of values.
156 * \param[in,out] values Values to convert.
157 * \param[in] type Type to convert to.
158 * \param[in] scanner Scanner data structure.
159 * \returns 0 on success, a non-zero value on error.
162 convert_values(t_selexpr_value
*values
, e_selvalue_t type
, void *scanner
)
164 t_selexpr_value
*value
;
171 rc1
= convert_value(value
, type
, scanner
);
172 if (rc1
!= 0 && rc
== 0)
178 /* FIXME: More informative error messages */
183 * Adds a child element for a parameter, keeping the parameter order.
185 * \param[in,out] root Root element to which the child is added.
186 * \param[in] child Child to add.
187 * \param[in] param Parameter for which this child is a value.
189 * Puts \p child in the child list of \p root such that the list remains
190 * in the same order as the corresponding parameters.
193 place_child(t_selelem
*root
, t_selelem
*child
, gmx_ana_selparam_t
*param
)
195 gmx_ana_selparam_t
*ps
;
198 ps
= root
->u
.expr
.method
->param
;
200 /* Put the child element in the correct place */
201 if (!root
->child
|| n
< root
->child
->u
.param
- ps
)
203 child
->next
= root
->child
;
211 while (prev
->next
&& prev
->next
->u
.param
- ps
>= n
)
215 child
->next
= prev
->next
;
221 * Comparison function for sorting integer ranges.
223 * \param[in] a Pointer to the first range.
224 * \param[in] b Pointer to the second range.
225 * \returns -1, 0, or 1 depending on the relative order of \p a and \p b.
227 * The ranges are primarily sorted based on their starting point, and
228 * secondarily based on length (longer ranges come first).
231 cmp_int_range(const void *a
, const void *b
)
233 if (((int *)a
)[0] < ((int *)b
)[0])
237 if (((int *)a
)[0] > ((int *)b
)[0])
241 if (((int *)a
)[1] > ((int *)b
)[1])
249 * Comparison function for sorting real ranges.
251 * \param[in] a Pointer to the first range.
252 * \param[in] b Pointer to the second range.
253 * \returns -1, 0, or 1 depending on the relative order of \p a and \p b.
255 * The ranges are primarily sorted based on their starting point, and
256 * secondarily based on length (longer ranges come first).
259 cmp_real_range(const void *a
, const void *b
)
261 if (((real
*)a
)[0] < ((real
*)b
)[0])
265 if (((real
*)a
)[0] > ((real
*)b
)[0])
269 if (((real
*)a
)[1] > ((real
*)b
)[1])
277 * Parses the values for a parameter that takes integer or real ranges.
279 * \param[in] nval Number of values in \p values.
280 * \param[in] values Pointer to the list of values.
281 * \param param Parameter to parse.
282 * \returns TRUE if the values were parsed successfully, FALSE otherwise.
285 parse_values_range(int nval
, t_selexpr_value
*values
, gmx_ana_selparam_t
*param
)
287 t_selexpr_value
*value
;
292 param
->flags
&= ~SPAR_DYNAMIC
;
293 if (param
->val
.type
!= INT_VALUE
&& param
->val
.type
!= REAL_VALUE
)
295 gmx_bug("internal error");
300 if (param
->val
.type
== INT_VALUE
)
314 _gmx_selparser_error("expressions not supported within range parameters");
317 if (value
->type
!= param
->val
.type
)
319 gmx_bug("internal error");
322 if (param
->val
.type
== INT_VALUE
)
324 /* Make sure the input range is in increasing order */
325 if (value
->u
.i
.i1
> value
->u
.i
.i2
)
327 int tmp
= value
->u
.i
.i1
;
328 value
->u
.i
.i1
= value
->u
.i
.i2
;
331 /* Check if the new range overlaps or extends the previous one */
332 if (i
> 0 && value
->u
.i
.i1
<= idata
[i
-1]+1 && value
->u
.i
.i2
>= idata
[i
-2]-1)
334 idata
[i
-2] = min(idata
[i
-2], value
->u
.i
.i1
);
335 idata
[i
-1] = max(idata
[i
-1], value
->u
.i
.i2
);
339 idata
[i
++] = value
->u
.i
.i1
;
340 idata
[i
++] = value
->u
.i
.i2
;
345 /* Make sure the input range is in increasing order */
346 if (value
->u
.r
.r1
> value
->u
.r
.r2
)
348 real tmp
= value
->u
.r
.r1
;
349 value
->u
.r
.r1
= value
->u
.r
.r2
;
352 /* Check if the new range overlaps or extends the previous one */
353 if (i
> 0 && value
->u
.r
.r1
<= rdata
[i
-1] && value
->u
.r
.r2
>= rdata
[i
-2])
355 rdata
[i
-2] = min(rdata
[i
-2], value
->u
.r
.r1
);
356 rdata
[i
-1] = max(rdata
[i
-1], value
->u
.r
.r2
);
360 rdata
[i
++] = value
->u
.r
.r1
;
361 rdata
[i
++] = value
->u
.r
.r2
;
367 /* Sort the ranges and merge consequent ones */
368 if (param
->val
.type
== INT_VALUE
)
370 qsort(idata
, n
, 2*sizeof(int), &cmp_int_range
);
371 for (i
= j
= 2; i
< 2*n
; i
+= 2)
373 if (idata
[j
-1]+1 >= idata
[i
])
375 if (idata
[i
+1] > idata
[j
-1])
377 idata
[j
-1] = idata
[i
+1];
383 idata
[j
+1] = idata
[i
+1];
390 qsort(rdata
, n
, 2*sizeof(real
), &cmp_real_range
);
391 for (i
= j
= 2; i
< 2*n
; i
+= 2)
393 if (rdata
[j
-1]+1 >= rdata
[i
])
395 if (rdata
[i
+1] > rdata
[j
-1])
397 rdata
[j
-1] = rdata
[i
+1];
403 rdata
[j
+1] = rdata
[i
+1];
409 /* Store the values */
410 if (param
->flags
& SPAR_VARNUM
)
413 if (param
->val
.type
== INT_VALUE
)
416 _gmx_selvalue_setstore_alloc(¶m
->val
, idata
, j
);
421 _gmx_selvalue_setstore_alloc(¶m
->val
, rdata
, j
);
426 if (n
!= param
->val
.nr
)
428 _gmx_selparser_error("the value of parameter '%s' should consist of exactly one range",
434 if (param
->val
.type
== INT_VALUE
)
436 memcpy(param
->val
.u
.i
, idata
, 2*n
*sizeof(int));
441 memcpy(param
->val
.u
.r
, rdata
, 2*n
*sizeof(real
));
447 *param
->nvalptr
= param
->val
.nr
;
449 param
->nvalptr
= NULL
;
455 * Parses the values for a parameter that takes a variable number of values.
457 * \param[in] nval Number of values in \p values.
458 * \param[in] values Pointer to the list of values.
459 * \param param Parameter to parse.
460 * \param root Selection element to which child expressions are added.
461 * \returns TRUE if the values were parsed successfully, FALSE otherwise.
463 * For integer ranges, the sequence of numbers from the first to second value
464 * is stored, each as a separate value.
467 parse_values_varnum(int nval
, t_selexpr_value
*values
,
468 gmx_ana_selparam_t
*param
, t_selelem
*root
)
470 t_selexpr_value
*value
;
473 param
->flags
&= ~SPAR_DYNAMIC
;
474 /* Update nval if there are integer ranges. */
475 if (param
->val
.type
== INT_VALUE
)
480 if (value
->type
== INT_VALUE
&& !value
->bExpr
)
482 nval
+= abs(value
->u
.i
.i2
- value
->u
.i
.i1
);
488 /* Check that the value type is actually implemented */
489 if (param
->val
.type
!= INT_VALUE
&& param
->val
.type
!= REAL_VALUE
490 && param
->val
.type
!= STR_VALUE
&& param
->val
.type
!= POS_VALUE
)
492 gmx_bug("internal error");
496 /* Reserve appropriate amount of memory */
497 if (param
->val
.type
== POS_VALUE
)
499 gmx_ana_pos_reserve(param
->val
.u
.p
, nval
, 0);
500 gmx_ana_pos_set_nr(param
->val
.u
.p
, nval
);
501 gmx_ana_indexmap_init(¶m
->val
.u
.p
->m
, NULL
, NULL
, INDEX_UNKNOWN
);
505 _gmx_selvalue_reserve(¶m
->val
, nval
);
514 _gmx_selparser_error("expressions not supported within value lists");
517 if (value
->type
!= param
->val
.type
)
519 gmx_bug("internal error");
522 switch (param
->val
.type
)
525 if (value
->u
.i
.i1
<= value
->u
.i
.i2
)
527 for (j
= value
->u
.i
.i1
; j
<= value
->u
.i
.i2
; ++j
)
529 param
->val
.u
.i
[i
++] = j
;
534 for (j
= value
->u
.i
.i1
; j
>= value
->u
.i
.i2
; --j
)
536 param
->val
.u
.i
[i
++] = j
;
541 if (value
->u
.r
.r1
!= value
->u
.r
.r2
)
543 _gmx_selparser_error("real ranges not supported for parameter '%s'", param
->name
);
546 param
->val
.u
.r
[i
++] = value
->u
.r
.r1
;
548 case STR_VALUE
: param
->val
.u
.s
[i
++] = strdup(value
->u
.s
); break;
549 case POS_VALUE
: copy_rvec(value
->u
.x
, param
->val
.u
.p
->x
[i
++]); break;
550 default: /* Should not be reached */
551 gmx_bug("internal error");
559 *param
->nvalptr
= param
->val
.nr
;
561 param
->nvalptr
= NULL
;
562 /* Create a dummy child element to store the string values.
563 * This element is responsible for freeing the values, but carries no
565 if (param
->val
.type
== STR_VALUE
)
569 child
= _gmx_selelem_create(SEL_CONST
);
570 _gmx_selelem_set_vtype(child
, STR_VALUE
);
571 child
->name
= param
->name
;
572 child
->flags
&= ~SEL_ALLOCVAL
;
573 child
->flags
|= SEL_FLAGSSET
| SEL_VARNUMVAL
| SEL_ALLOCDATA
;
574 child
->v
.nr
= param
->val
.nr
;
575 _gmx_selvalue_setstore(&child
->v
, param
->val
.u
.s
);
576 /* Because the child is not group-valued, the u union is not used
577 * for anything, so we can abuse it by storing the parameter value
578 * as place_child() expects, but this is really ugly... */
579 child
->u
.param
= param
;
580 place_child(root
, child
, param
);
587 * Adds a new subexpression reference to a selection element.
589 * \param[in,out] root Root element to which the subexpression is added.
590 * \param[in] param Parameter for which this expression is a value.
591 * \param[in] expr Expression to add.
592 * \returns The created child element.
594 * Creates a new \ref SEL_SUBEXPRREF element and adds it into the child
596 * If \p expr is already a \ref SEL_SUBEXPRREF, it is used as it is.
597 * \ref SEL_ALLOCVAL is cleared for the returned element.
600 add_child(t_selelem
*root
, gmx_ana_selparam_t
*param
, t_selelem
*expr
)
605 if (root
->type
!= SEL_EXPRESSION
&& root
->type
!= SEL_MODIFIER
)
607 gmx_bug("unsupported root element for selection parameter parser");
610 /* Create a subexpression reference element if necessary */
611 if (expr
->type
== SEL_SUBEXPRREF
)
617 child
= _gmx_selelem_create(SEL_SUBEXPRREF
);
622 _gmx_selelem_set_vtype(child
, expr
->v
.type
);
625 /* Setup the child element */
626 child
->flags
&= ~SEL_ALLOCVAL
;
627 child
->u
.param
= param
;
628 if (child
->v
.type
!= param
->val
.type
)
630 _gmx_selparser_error("invalid expression value for parameter '%s'",
634 rc
= _gmx_selelem_update_flags(child
);
639 if ((child
->flags
& SEL_DYNAMIC
) && !(param
->flags
& SPAR_DYNAMIC
))
641 _gmx_selparser_error("parameter '%s' does not support dynamic values",
645 if (!(child
->flags
& SEL_DYNAMIC
))
647 param
->flags
&= ~SPAR_DYNAMIC
;
649 /* Put the child element in the correct place */
650 place_child(root
, child
, param
);
656 _gmx_selelem_free(child
);
662 * Parses an expression value for a parameter that takes a variable number of values.
664 * \param[in] nval Number of values in \p values.
665 * \param[in] values Pointer to the list of values.
666 * \param param Parameter to parse.
667 * \param root Selection element to which child expressions are added.
668 * \returns TRUE if the values were parsed successfully, FALSE otherwise.
671 parse_values_varnum_expr(int nval
, t_selexpr_value
*values
,
672 gmx_ana_selparam_t
*param
, t_selelem
*root
)
674 t_selexpr_value
*value
;
678 if (nval
!= 1 || !values
->bExpr
)
680 gmx_bug("internal error");
685 child
= add_child(root
, param
, value
->u
.expr
);
686 value
->u
.expr
= NULL
;
692 /* Process single-valued expressions */
693 /* TODO: We should also handle SEL_SINGLEVAL expressions here */
694 if (child
->v
.type
== POS_VALUE
|| child
->v
.type
== GROUP_VALUE
)
696 /* Set the value storage */
697 _gmx_selvalue_setstore(&child
->v
, param
->val
.u
.ptr
);
701 *param
->nvalptr
= param
->val
.nr
;
703 param
->nvalptr
= NULL
;
707 if (!(child
->flags
& SEL_VARNUMVAL
))
709 _gmx_selparser_error("invalid expression value for parameter '%s'",
714 child
->flags
|= SEL_ALLOCVAL
;
716 *param
->nvalptr
= param
->val
.nr
;
717 /* Rest of the initialization is done during compilation in
724 * Initializes the storage of an expression value.
726 * \param[in,out] sel Selection element that evaluates the value.
727 * \param[in] param Parameter to receive the value.
728 * \param[in] i The value of \p sel evaluates the value \p i for
731 * Initializes the data pointer of \p sel such that the result is stored
732 * as the value \p i of \p param.
733 * This function is used internally by parse_values_std().
736 set_expr_value_store(t_selelem
*sel
, gmx_ana_selparam_t
*param
, int i
)
738 if (sel
->v
.type
!= GROUP_VALUE
&& !(sel
->flags
& SEL_SINGLEVAL
))
740 _gmx_selparser_error("invalid expression value for parameter '%s'",
746 case INT_VALUE
: sel
->v
.u
.i
= ¶m
->val
.u
.i
[i
]; break;
747 case REAL_VALUE
: sel
->v
.u
.r
= ¶m
->val
.u
.r
[i
]; break;
748 case STR_VALUE
: sel
->v
.u
.s
= ¶m
->val
.u
.s
[i
]; break;
749 case POS_VALUE
: sel
->v
.u
.p
= ¶m
->val
.u
.p
[i
]; break;
750 case GROUP_VALUE
: sel
->v
.u
.g
= ¶m
->val
.u
.g
[i
]; break;
752 gmx_bug("internal error");
761 * Parses the values for a parameter that takes a constant number of values.
763 * \param[in] nval Number of values in \p values.
764 * \param[in] values Pointer to the list of values.
765 * \param param Parameter to parse.
766 * \param root Selection element to which child expressions are added.
767 * \returns TRUE if the values were parsed successfully, FALSE otherwise.
769 * For integer ranges, the sequence of numbers from the first to second value
770 * is stored, each as a separate value.
773 parse_values_std(int nval
, t_selexpr_value
*values
, gmx_ana_selparam_t
*param
,
776 t_selexpr_value
*value
;
781 /* Handle atom-valued parameters */
782 if (param
->flags
& SPAR_ATOMVAL
)
786 _gmx_selparser_error("extra values for parameter '%s'", param
->name
);
792 child
= add_child(root
, param
, value
->u
.expr
);
793 value
->u
.expr
= NULL
;
798 child
->flags
|= SEL_ALLOCVAL
;
799 if (child
->v
.type
!= GROUP_VALUE
&& (child
->flags
& SEL_ATOMVAL
))
801 /* Rest of the initialization is done during compilation in
803 /* TODO: Positions are not correctly handled */
807 *param
->nvalptr
= -1;
811 param
->flags
&= ~SPAR_ATOMVAL
;
817 param
->nvalptr
= NULL
;
818 if (param
->val
.type
== INT_VALUE
|| param
->val
.type
== REAL_VALUE
819 || param
->val
.type
== STR_VALUE
)
821 _gmx_selvalue_reserve(¶m
->val
, 1);
823 return set_expr_value_store(child
, param
, 0);
825 /* If we reach here, proceed with normal parameter handling */
827 if (param
->val
.type
== INT_VALUE
|| param
->val
.type
== REAL_VALUE
828 || param
->val
.type
== STR_VALUE
)
830 _gmx_selvalue_reserve(¶m
->val
, 1);
832 param
->flags
&= ~SPAR_ATOMVAL
;
833 param
->flags
&= ~SPAR_DYNAMIC
;
839 while (value
&& i
< param
->val
.nr
)
841 if (value
->type
!= param
->val
.type
)
843 _gmx_selparser_error("incorrect value for parameter '%s' skipped", param
->name
);
849 child
= add_child(root
, param
, value
->u
.expr
);
850 /* Clear the expression from the value once it is stored */
851 value
->u
.expr
= NULL
;
852 /* Check that the expression is valid */
857 if (!set_expr_value_store(child
, param
, i
))
861 if (child
->flags
& SEL_DYNAMIC
)
868 /* Value is not an expression */
872 if (value
->u
.i
.i1
<= value
->u
.i
.i2
)
874 for (j
= value
->u
.i
.i1
; j
<= value
->u
.i
.i2
&& i
< param
->val
.nr
; ++j
)
876 param
->val
.u
.i
[i
++] = j
;
878 if (j
!= value
->u
.i
.i2
+ 1)
880 _gmx_selparser_error("extra values for parameter '%s' skipped", param
->name
);
885 for (j
= value
->u
.i
.i1
; j
>= value
->u
.i
.i2
&& i
< param
->val
.nr
; --j
)
887 param
->val
.u
.i
[i
++] = j
;
889 if (j
!= value
->u
.i
.i2
- 1)
891 _gmx_selparser_error("extra values for parameter '%s' skipped", param
->name
);
897 if (value
->u
.r
.r1
!= value
->u
.r
.r2
)
899 _gmx_selparser_error("real ranges not supported for parameter '%s'", param
->name
);
902 param
->val
.u
.r
[i
] = value
->u
.r
.r1
;
905 param
->val
.u
.s
[i
] = strdup(value
->u
.s
);
908 gmx_ana_pos_init_const(¶m
->val
.u
.p
[i
], value
->u
.x
);
912 gmx_bug("internal error");
921 _gmx_selparser_error("extra values for parameter '%s'", param
->name
);
924 if (i
< param
->val
.nr
)
926 _gmx_selparser_error("not enough values for parameter '%s'", param
->name
);
931 param
->flags
&= ~SPAR_DYNAMIC
;
935 *param
->nvalptr
= param
->val
.nr
;
937 param
->nvalptr
= NULL
;
943 * Parses the values for a gmx_boolean parameter.
945 * \param[in] name Name by which the parameter was given.
946 * \param[in] nval Number of values in \p values.
947 * \param[in] values Pointer to the list of values.
948 * \param param Parameter to parse.
949 * \returns TRUE if the values were parsed successfully, FALSE otherwise.
952 parse_values_gmx_bool(const char *name
, int nval
, t_selexpr_value
*values
, gmx_ana_selparam_t
*param
)
957 if (param
->val
.type
!= NO_VALUE
)
959 gmx_bug("internal error");
962 if (nval
> 1 || (values
&& values
->type
!= INT_VALUE
))
964 _gmx_selparser_error("gmx_boolean parameter '%s' takes only a yes/no/on/off/0/1 value", param
->name
);
969 /* Check if the parameter name is given with a 'no' prefix */
971 if (len
> 2 && name
[0] == 'n' && name
[1] == 'o'
972 && strncmp(name
+2, param
->name
, len
-2) == 0)
976 if (bSetNo
&& nval
> 0)
978 _gmx_selparser_error("gmx_boolean parameter 'no%s' should not have a value", param
->name
);
981 if (values
&& values
->u
.i
.i1
== 0)
986 *param
->val
.u
.b
= bSetNo
? FALSE
: TRUE
;
991 * Parses the values for an enumeration parameter.
993 * \param[in] nval Number of values in \p values.
994 * \param[in] values Pointer to the list of values.
995 * \param param Parameter to parse.
996 * \returns TRUE if the values were parsed successfully, FALSE otherwise.
999 parse_values_enum(int nval
, t_selexpr_value
*values
, gmx_ana_selparam_t
*param
)
1005 _gmx_selparser_error("a single value is required for parameter '%s'", param
->name
);
1008 if (values
->type
!= STR_VALUE
|| param
->val
.type
!= STR_VALUE
)
1010 gmx_bug("internal error");
1015 _gmx_selparser_error("expression value for enumerated parameter '%s' not supported", param
->name
);
1019 len
= strlen(values
->u
.s
);
1022 while (param
->val
.u
.s
[i
] != NULL
)
1024 if (strncmp(values
->u
.s
, param
->val
.u
.s
[i
], len
) == 0)
1026 /* Check if there is a duplicate match */
1029 _gmx_selparser_error("ambiguous value for parameter '%s'", param
->name
);
1038 _gmx_selparser_error("invalid value for parameter '%s'", param
->name
);
1041 param
->val
.u
.s
[0] = param
->val
.u
.s
[match
];
1046 * Replaces constant expressions with their values.
1048 * \param[in,out] values First element in the value list to process.
1051 convert_const_values(t_selexpr_value
*values
)
1053 t_selexpr_value
*val
;
1058 if (val
->bExpr
&& val
->u
.expr
->v
.type
!= GROUP_VALUE
&&
1059 val
->u
.expr
->type
== SEL_CONST
)
1061 t_selelem
*expr
= val
->u
.expr
;
1063 switch (expr
->v
.type
)
1066 val
->u
.i
.i1
= val
->u
.i
.i2
= expr
->v
.u
.i
[0];
1069 val
->u
.r
.r1
= val
->u
.r
.r2
= expr
->v
.u
.r
[0];
1072 val
->u
.s
= expr
->v
.u
.s
[0];
1075 copy_rvec(expr
->v
.u
.p
->x
[0], val
->u
.x
);
1078 gmx_bug("internal error");
1081 _gmx_selelem_free(expr
);
1088 * \param pparams List of parameters from the selection parser.
1089 * \param[in] nparam Number of parameters in \p params.
1090 * \param params Array of parameters to parse.
1091 * \param root Selection element to which child expressions are added.
1092 * \param[in] scanner Scanner data structure.
1093 * \returns TRUE if the parameters were parsed successfully, FALSE otherwise.
1095 * Initializes the \p params array based on the parameters in \p pparams.
1096 * See the documentation of \c gmx_ana_selparam_t for different options
1097 * available for parsing.
1099 * The list \p pparams and any associated values are freed after the parameters
1100 * have been processed, no matter is there was an error or not.
1103 _gmx_sel_parse_params(t_selexpr_param
*pparams
, int nparam
, gmx_ana_selparam_t
*params
,
1104 t_selelem
*root
, void *scanner
)
1106 t_selexpr_param
*pparam
;
1107 gmx_ana_selparam_t
*oparam
;
1111 /* Check that the value pointers of SPAR_VARNUM parameters are NULL and
1112 * that they are not NULL for other parameters */
1114 for (i
= 0; i
< nparam
; ++i
)
1116 if (params
[i
].val
.type
!= POS_VALUE
&& (params
[i
].flags
& (SPAR_VARNUM
| SPAR_ATOMVAL
)))
1118 if (params
[i
].val
.u
.ptr
!= NULL
)
1120 _gmx_selparser_error("warning: value pointer of parameter '%s' is not NULL\n"
1121 " although it should be for SPAR_VARNUM and SPAR_ATOMVAL parameters\n",
1124 if ((params
[i
].flags
& SPAR_VARNUM
)
1125 && (params
[i
].flags
& SPAR_DYNAMIC
) && !params
[i
].nvalptr
)
1127 _gmx_selparser_error("error: nvalptr of parameter '%s' is NULL\n"
1128 " but both SPAR_VARNUM and SPAR_DYNAMIC are specified\n",
1135 if (params
[i
].val
.u
.ptr
== NULL
)
1137 _gmx_selparser_error("error: value pointer of parameter '%s' is NULL\n",
1145 _gmx_selexpr_free_params(pparams
);
1148 /* Parse the parameters */
1153 /* Find the parameter and make some checks */
1154 if (pparam
->name
!= NULL
)
1157 oparam
= gmx_ana_selparam_find(pparam
->name
, nparam
, params
);
1161 oparam
= ¶ms
[i
];
1162 if (oparam
->name
!= NULL
)
1165 _gmx_selparser_error("too many NULL parameters provided");
1173 _gmx_selparser_error("all NULL parameters should appear in the beginning of the list");
1175 pparam
= pparam
->next
;
1180 _gmx_selparser_error("unknown parameter '%s' skipped", pparam
->name
);
1184 if (oparam
->flags
& SPAR_SET
)
1186 _gmx_selparser_error("parameter '%s' set multiple times, extra values skipped", pparam
->name
);
1190 oparam
->flags
|= SPAR_SET
;
1191 /* Process the values for the parameter */
1192 convert_const_values(pparam
->value
);
1193 if (convert_values(pparam
->value
, oparam
->val
.type
, scanner
) != 0)
1195 _gmx_selparser_error("invalid value for parameter '%s'", pparam
->name
);
1199 if (oparam
->val
.type
== NO_VALUE
)
1201 rc
= parse_values_gmx_bool(pparam
->name
, pparam
->nval
, pparam
->value
, oparam
);
1203 else if (oparam
->flags
& SPAR_RANGES
)
1205 rc
= parse_values_range(pparam
->nval
, pparam
->value
, oparam
);
1207 else if (oparam
->flags
& SPAR_VARNUM
)
1209 if (pparam
->nval
== 1 && pparam
->value
->bExpr
)
1211 rc
= parse_values_varnum_expr(pparam
->nval
, pparam
->value
, oparam
, root
);
1215 rc
= parse_values_varnum(pparam
->nval
, pparam
->value
, oparam
, root
);
1218 else if (oparam
->flags
& SPAR_ENUMVAL
)
1220 rc
= parse_values_enum(pparam
->nval
, pparam
->value
, oparam
);
1224 rc
= parse_values_std(pparam
->nval
, pparam
->value
, oparam
, root
);
1230 /* Advance to the next parameter */
1232 pparam
= pparam
->next
;
1234 /* Check that all required parameters are present */
1235 for (i
= 0; i
< nparam
; ++i
)
1237 if (!(params
[i
].flags
& SPAR_OPTIONAL
) && !(params
[i
].flags
& SPAR_SET
))
1239 _gmx_selparser_error("required parameter '%s' not specified", params
[i
].name
);
1244 _gmx_selexpr_free_params(pparams
);