Remove some hardcoded selection string buffers
[gromacs.git] / src / gromacs / selection / params.cpp
blobe0305da77297c723957ef572e65a65e5ae269e96
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2009,2010,2011,2012,2013,2014, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
35 /*! \internal \file
36 * \brief
37 * Implements functions in selparam.h.
39 * \author Teemu Murtola <teemu.murtola@gmail.com>
40 * \ingroup module_selection
42 #include <algorithm>
43 #include <string>
45 #include "gromacs/legacyheaders/vec.h"
47 #include "gromacs/selection/position.h"
48 #include "gromacs/selection/selmethod.h"
49 #include "gromacs/selection/selparam.h"
50 #include "gromacs/utility/exceptions.h"
51 #include "gromacs/utility/gmxassert.h"
52 #include "gromacs/utility/messagestringcollector.h"
53 #include "gromacs/utility/smalloc.h"
54 #include "gromacs/utility/stringutil.h"
56 #include "parsetree.h"
57 #include "position.h"
58 #include "scanner.h"
59 #include "selelem.h"
61 using gmx::SelectionParserValue;
62 using gmx::SelectionParserValueList;
63 using gmx::SelectionParserParameter;
64 using gmx::SelectionParserParameterList;
65 using gmx::SelectionTreeElement;
66 using gmx::SelectionTreeElementPointer;
68 /*!
69 * \param[in] name Name of the parameter to search.
70 * \param[in] nparam Number of parameters in the \p param array.
71 * \param[in] param Parameter array to search.
72 * \returns Pointer to the parameter in the \p param
73 * or NULL if no parameter with name \p name was found.
75 * The comparison is case-sensitive.
77 gmx_ana_selparam_t *
78 gmx_ana_selparam_find(const char *name, int nparam, gmx_ana_selparam_t *param)
80 int i;
82 if (nparam == 0)
84 return NULL;
86 /* Find the first non-null parameter */
87 i = 0;
88 while (i < nparam && param[i].name == NULL)
90 ++i;
92 /* Process the special case of a NULL parameter */
93 if (name == NULL)
95 return (i == 0) ? NULL : &param[i-1];
97 for (; i < nparam; ++i)
99 if (!strcmp(param[i].name, name))
101 return &param[i];
103 /* Check for 'no' prefix on boolean parameters */
104 if (param[i].val.type == NO_VALUE
105 && strlen(name) > 2 && name[0] == 'n' && name[1] == 'o'
106 && !strcmp(param[i].name, name+2))
108 return &param[i];
111 return NULL;
114 /*! \brief
115 * Does a type conversion on a SelectionParserValue.
117 * \param[in,out] value Value to convert.
118 * \param[in] type Type to convert to.
119 * \param[in] scanner Scanner data structure.
120 * \returns 0 on success, a non-zero value on error.
122 static int
123 convert_value(SelectionParserValue *value, e_selvalue_t type, void *scanner)
125 if (value->type == type || type == NO_VALUE)
127 return 0;
129 if (value->hasExpressionValue())
131 /* Conversion from atom selection to position using default
132 * reference positions. */
133 if (value->type == GROUP_VALUE && type == POS_VALUE)
135 SelectionTreeElementPointer expr =
136 _gmx_sel_init_position(value->expr, NULL, scanner);
137 // FIXME: Use exceptions
138 if (!expr)
140 return -1;
142 *value = SelectionParserValue::createExpr(expr);
143 return 0;
145 return -1;
147 else
149 /* Integers to floating point are easy */
150 if (value->type == INT_VALUE && type == REAL_VALUE)
152 *value = SelectionParserValue::createRealRange(value->u.i.i1,
153 value->u.i.i2);
154 return 0;
156 /* Reals that are integer-valued can also be converted */
157 if (value->type == REAL_VALUE && type == INT_VALUE)
159 int i1 = static_cast<int>(value->u.r.r1);
160 int i2 = static_cast<int>(value->u.r.r2);
161 if (gmx_within_tol(value->u.r.r1, i1, GMX_REAL_EPS)
162 && gmx_within_tol(value->u.r.r2, i2, GMX_REAL_EPS))
164 *value = SelectionParserValue::createIntegerRange(i1, i2);
165 return 0;
169 return -1;
172 /*! \brief
173 * Does a type conversion on a list of values.
175 * \param[in,out] values Values to convert.
176 * \param[in] type Type to convert to.
177 * \param[in] scanner Scanner data structure.
178 * \returns 0 on success, a non-zero value on error.
180 static int
181 convert_values(SelectionParserValueList *values, e_selvalue_t type, void *scanner)
183 int rc = 0;
184 SelectionParserValueList::iterator value;
185 for (value = values->begin(); value != values->end(); ++value)
187 int rc1 = convert_value(&*value, type, scanner);
188 if (rc1 != 0 && rc == 0)
190 rc = rc1;
193 /* FIXME: More informative error messages */
194 return rc;
197 /*! \brief
198 * Adds a child element for a parameter, keeping the parameter order.
200 * \param[in,out] root Root element to which the child is added.
201 * \param[in] child Child to add.
202 * \param[in] param Parameter for which this child is a value.
204 * Puts \p child in the child list of \p root such that the list remains
205 * in the same order as the corresponding parameters.
207 static void
208 place_child(const SelectionTreeElementPointer &root,
209 const SelectionTreeElementPointer &child,
210 gmx_ana_selparam_t *param)
212 gmx_ana_selparam_t *ps;
213 int n;
215 ps = root->u.expr.method->param;
216 n = param - ps;
217 /* Put the child element in the correct place */
218 if (!root->child || n < root->child->u.param - ps)
220 child->next = root->child;
221 root->child = child;
223 else
225 SelectionTreeElementPointer prev = root->child;
226 while (prev->next && prev->next->u.param - ps >= n)
228 prev = prev->next;
230 child->next = prev->next;
231 prev->next = child;
235 /*! \brief
236 * Comparison function for sorting integer ranges.
238 * \param[in] a Pointer to the first range.
239 * \param[in] b Pointer to the second range.
240 * \returns -1, 0, or 1 depending on the relative order of \p a and \p b.
242 * The ranges are primarily sorted based on their starting point, and
243 * secondarily based on length (longer ranges come first).
245 static int
246 cmp_int_range(const void *a, const void *b)
248 if (((int *)a)[0] < ((int *)b)[0])
250 return -1;
252 if (((int *)a)[0] > ((int *)b)[0])
254 return 1;
256 if (((int *)a)[1] > ((int *)b)[1])
258 return -1;
260 return 0;
263 /*! \brief
264 * Comparison function for sorting real ranges.
266 * \param[in] a Pointer to the first range.
267 * \param[in] b Pointer to the second range.
268 * \returns -1, 0, or 1 depending on the relative order of \p a and \p b.
270 * The ranges are primarily sorted based on their starting point, and
271 * secondarily based on length (longer ranges come first).
273 static int
274 cmp_real_range(const void *a, const void *b)
276 if (((real *)a)[0] < ((real *)b)[0])
278 return -1;
280 if (((real *)a)[0] > ((real *)b)[0])
282 return 1;
284 if (((real *)a)[1] > ((real *)b)[1])
286 return -1;
288 return 0;
291 /*! \brief
292 * Parses the values for a parameter that takes integer or real ranges.
294 * \param[in] values List of values.
295 * \param param Parameter to parse.
296 * \param[in] scanner Scanner data structure.
297 * \returns true if the values were parsed successfully, false otherwise.
299 static bool
300 parse_values_range(const SelectionParserValueList &values,
301 gmx_ana_selparam_t *param, void *scanner)
303 int *idata;
304 real *rdata;
305 int i, j, n;
307 param->flags &= ~SPAR_DYNAMIC;
308 GMX_RELEASE_ASSERT(param->val.type == INT_VALUE || param->val.type == REAL_VALUE,
309 "Invalid range parameter type");
310 idata = NULL;
311 rdata = NULL;
312 if (param->val.type == INT_VALUE)
314 snew(idata, values.size()*2);
316 else
318 snew(rdata, values.size()*2);
320 i = 0;
321 SelectionParserValueList::const_iterator value;
322 for (value = values.begin(); value != values.end(); ++value)
324 if (value->hasExpressionValue())
326 _gmx_selparser_error(scanner, "expressions not supported within range parameters");
327 return false;
329 GMX_RELEASE_ASSERT(value->type == param->val.type,
330 "Invalid range value type (should have been caught earlier)");
331 if (param->val.type == INT_VALUE)
333 int i1 = std::min(value->u.i.i1, value->u.i.i2);
334 int i2 = std::max(value->u.i.i1, value->u.i.i2);
335 /* Check if the new range overlaps or extends the previous one */
336 if (i > 0 && i1 <= idata[i-1]+1 && i2 >= idata[i-2]-1)
338 idata[i-2] = std::min(idata[i-2], i1);
339 idata[i-1] = std::max(idata[i-1], i2);
341 else
343 idata[i++] = i1;
344 idata[i++] = i2;
347 else
349 real r1 = std::min(value->u.r.r1, value->u.r.r2);
350 real r2 = std::max(value->u.r.r1, value->u.r.r2);
351 /* Check if the new range overlaps or extends the previous one */
352 if (i > 0 && r1 <= rdata[i-1] && r2 >= rdata[i-2])
354 rdata[i-2] = std::min(rdata[i-2], r1);
355 rdata[i-1] = std::max(rdata[i-1], r2);
357 else
359 rdata[i++] = r1;
360 rdata[i++] = r2;
364 n = i/2;
365 /* Sort the ranges and merge consequent ones */
366 if (param->val.type == INT_VALUE)
368 qsort(idata, n, 2*sizeof(int), &cmp_int_range);
369 for (i = j = 2; i < 2*n; i += 2)
371 if (idata[j-1]+1 >= idata[i])
373 if (idata[i+1] > idata[j-1])
375 idata[j-1] = idata[i+1];
378 else
380 idata[j] = idata[i];
381 idata[j+1] = idata[i+1];
382 j += 2;
386 else
388 qsort(rdata, n, 2*sizeof(real), &cmp_real_range);
389 for (i = j = 2; i < 2*n; i += 2)
391 if (rdata[j-1]+1 >= rdata[i])
393 if (rdata[i+1] > rdata[j-1])
395 rdata[j-1] = rdata[i+1];
398 else
400 rdata[j] = rdata[i];
401 rdata[j+1] = rdata[i+1];
402 j += 2;
406 n = j/2;
407 /* Store the values */
408 if (param->flags & SPAR_VARNUM)
410 param->val.nr = n;
411 if (param->val.type == INT_VALUE)
413 srenew(idata, j);
414 _gmx_selvalue_setstore_alloc(&param->val, idata, j);
416 else
418 srenew(rdata, j);
419 _gmx_selvalue_setstore_alloc(&param->val, rdata, j);
422 else
424 if (n != param->val.nr)
426 _gmx_selparser_error(scanner, "the value should consist of exactly one range");
427 sfree(idata);
428 sfree(rdata);
429 return false;
431 if (param->val.type == INT_VALUE)
433 memcpy(param->val.u.i, idata, 2*n*sizeof(int));
434 sfree(idata);
436 else
438 memcpy(param->val.u.r, rdata, 2*n*sizeof(real));
439 sfree(rdata);
442 if (param->nvalptr)
444 *param->nvalptr = param->val.nr;
446 param->nvalptr = NULL;
448 return true;
451 /*! \brief
452 * Parses the values for a parameter that takes a variable number of values.
454 * \param[in] values List of values.
455 * \param param Parameter to parse.
456 * \param root Selection element to which child expressions are added.
457 * \param[in] scanner Scanner data structure.
458 * \returns true if the values were parsed successfully, false otherwise.
460 * For integer ranges, the sequence of numbers from the first to second value
461 * is stored, each as a separate value.
463 static bool
464 parse_values_varnum(const SelectionParserValueList &values,
465 gmx_ana_selparam_t *param,
466 const SelectionTreeElementPointer &root,
467 void *scanner)
469 int i, j;
471 param->flags &= ~SPAR_DYNAMIC;
472 /* Compute number of values, considering also integer ranges. */
473 size_t valueCount = values.size();
474 if (param->val.type == INT_VALUE)
476 SelectionParserValueList::const_iterator value;
477 for (value = values.begin(); value != values.end(); ++value)
479 if (value->type == INT_VALUE && !value->hasExpressionValue())
481 valueCount += abs(value->u.i.i2 - value->u.i.i1);
486 /* Check that the value type is actually implemented */
487 if (param->val.type != INT_VALUE && param->val.type != REAL_VALUE
488 && param->val.type != STR_VALUE && param->val.type != POS_VALUE)
490 GMX_THROW(gmx::InternalError("Variable-count value type not implemented"));
493 /* Reserve appropriate amount of memory */
494 if (param->val.type == POS_VALUE)
496 gmx_ana_pos_reserve(param->val.u.p, valueCount, 0);
497 gmx_ana_indexmap_init(&param->val.u.p->m, NULL, NULL, INDEX_UNKNOWN);
498 gmx_ana_pos_set_nr(param->val.u.p, valueCount);
500 else
502 _gmx_selvalue_reserve(&param->val, valueCount);
505 i = 0;
506 SelectionParserValueList::const_iterator value;
507 for (value = values.begin(); value != values.end(); ++value)
509 if (value->hasExpressionValue())
511 _gmx_selparser_error(scanner, "expressions not supported within value lists");
512 return false;
514 GMX_RELEASE_ASSERT(value->type == param->val.type,
515 "Invalid value type (should have been caught earlier)");
516 switch (param->val.type)
518 case INT_VALUE:
519 if (value->u.i.i1 <= value->u.i.i2)
521 for (j = value->u.i.i1; j <= value->u.i.i2; ++j)
523 param->val.u.i[i++] = j;
526 else
528 for (j = value->u.i.i1; j >= value->u.i.i2; --j)
530 param->val.u.i[i++] = j;
533 break;
534 case REAL_VALUE:
535 if (value->u.r.r1 != value->u.r.r2)
537 _gmx_selparser_error(scanner, "real ranges not supported");
538 return false;
540 param->val.u.r[i++] = value->u.r.r1;
541 break;
542 case STR_VALUE:
543 param->val.u.s[i++] = strdup(value->stringValue().c_str());
544 break;
545 case POS_VALUE: copy_rvec(value->u.x, param->val.u.p->x[i++]); break;
546 default: /* Should not be reached */
547 GMX_THROW(gmx::InternalError("Variable-count value type not implemented"));
550 param->val.nr = i;
551 if (param->nvalptr)
553 *param->nvalptr = param->val.nr;
555 param->nvalptr = NULL;
556 /* Create a dummy child element to store the string values.
557 * This element is responsible for freeing the values, but carries no
558 * other function. */
559 if (param->val.type == STR_VALUE)
561 SelectionTreeElementPointer child(new SelectionTreeElement(SEL_CONST));
562 _gmx_selelem_set_vtype(child, STR_VALUE);
563 child->setName(param->name);
564 child->flags &= ~SEL_ALLOCVAL;
565 child->flags |= SEL_FLAGSSET | SEL_VARNUMVAL | SEL_ALLOCDATA;
566 child->v.nr = param->val.nr;
567 _gmx_selvalue_setstore(&child->v, param->val.u.s);
568 /* Because the child is not group-valued, the u union is not used
569 * for anything, so we can abuse it by storing the parameter value
570 * as place_child() expects, but this is really ugly... */
571 child->u.param = param;
572 place_child(root, child, param);
575 return true;
578 /*! \brief
579 * Adds a new subexpression reference to a selection element.
581 * \param[in,out] root Root element to which the subexpression is added.
582 * \param[in] param Parameter for which this expression is a value.
583 * \param[in] expr Expression to add.
584 * \param[in] scanner Scanner data structure.
585 * \returns The created child element.
587 * Creates a new \ref SEL_SUBEXPRREF element and adds it into the child
588 * list of \p root.
589 * If \p expr is already a \ref SEL_SUBEXPRREF, it is used as it is.
590 * \ref SEL_ALLOCVAL is cleared for the returned element.
592 static SelectionTreeElementPointer
593 add_child(const SelectionTreeElementPointer &root, gmx_ana_selparam_t *param,
594 const SelectionTreeElementPointer &expr, void *scanner)
596 GMX_RELEASE_ASSERT(root->type == SEL_EXPRESSION || root->type == SEL_MODIFIER,
597 "Unsupported root element for selection parameter parser");
598 SelectionTreeElementPointer child;
599 /* Create a subexpression reference element if necessary */
600 if (expr->type == SEL_SUBEXPRREF)
602 child = expr;
604 else
606 child.reset(new SelectionTreeElement(SEL_SUBEXPRREF));
607 _gmx_selelem_set_vtype(child, expr->v.type);
608 child->child = expr;
610 /* Setup the child element */
611 child->flags &= ~SEL_ALLOCVAL;
612 child->u.param = param;
613 if (child->v.type != param->val.type)
615 _gmx_selparser_error(scanner, "invalid expression value");
616 // FIXME: Use exceptions.
617 return SelectionTreeElementPointer();
619 _gmx_selelem_update_flags(child);
620 if ((child->flags & SEL_DYNAMIC) && !(param->flags & SPAR_DYNAMIC))
622 _gmx_selparser_error(scanner, "dynamic values not supported");
623 // FIXME: Use exceptions.
624 return SelectionTreeElementPointer();
626 if (!(child->flags & SEL_DYNAMIC))
628 param->flags &= ~SPAR_DYNAMIC;
630 /* Put the child element in the correct place */
631 place_child(root, child, param);
632 return child;
635 /*! \brief
636 * Parses an expression value for a parameter that takes a variable number of values.
638 * \param[in] values List of values.
639 * \param param Parameter to parse.
640 * \param root Selection element to which child expressions are added.
641 * \param[in] scanner Scanner data structure.
642 * \returns true if the values were parsed successfully, false otherwise.
644 static bool
645 parse_values_varnum_expr(const SelectionParserValueList &values,
646 gmx_ana_selparam_t *param,
647 const SelectionTreeElementPointer &root,
648 void *scanner)
650 GMX_RELEASE_ASSERT(values.size() == 1 && values.front().hasExpressionValue(),
651 "Called with an invalid type of value");
653 SelectionTreeElementPointer child
654 = add_child(root, param, values.front().expr, scanner);
655 if (!child)
657 return false;
660 /* Process single-valued expressions */
661 /* TODO: We should also handle SEL_SINGLEVAL expressions here */
662 if (child->v.type == POS_VALUE || child->v.type == GROUP_VALUE)
664 /* Set the value storage */
665 _gmx_selvalue_setstore(&child->v, param->val.u.ptr);
666 param->val.nr = 1;
667 if (param->nvalptr)
669 *param->nvalptr = param->val.nr;
671 param->nvalptr = NULL;
672 return true;
675 if (!(child->flags & SEL_VARNUMVAL))
677 _gmx_selparser_error(scanner, "invalid expression value");
678 return false;
681 child->flags |= SEL_ALLOCVAL;
682 param->val.nr = -1;
683 *param->nvalptr = param->val.nr;
684 /* Rest of the initialization is done during compilation in
685 * init_method(). */
687 return true;
690 /*! \brief
691 * Initializes the storage of an expression value.
693 * \param[in,out] sel Selection element that evaluates the value.
694 * \param[in] param Parameter to receive the value.
695 * \param[in] i The value of \p sel evaluates the value \p i for
696 * \p param.
697 * \param[in] scanner Scanner data structure.
699 * Initializes the data pointer of \p sel such that the result is stored
700 * as the value \p i of \p param.
701 * This function is used internally by parse_values_std().
703 static bool
704 set_expr_value_store(const SelectionTreeElementPointer &sel,
705 gmx_ana_selparam_t *param, int i, void *scanner)
707 if (sel->v.type != GROUP_VALUE && !(sel->flags & SEL_SINGLEVAL))
709 _gmx_selparser_error(scanner, "invalid expression value");
710 return false;
712 switch (sel->v.type)
714 case INT_VALUE: sel->v.u.i = &param->val.u.i[i]; break;
715 case REAL_VALUE: sel->v.u.r = &param->val.u.r[i]; break;
716 case STR_VALUE: sel->v.u.s = &param->val.u.s[i]; break;
717 case POS_VALUE: sel->v.u.p = &param->val.u.p[i]; break;
718 case GROUP_VALUE: sel->v.u.g = &param->val.u.g[i]; break;
719 default: /* Error */
720 GMX_THROW(gmx::InternalError("Invalid value type"));
722 sel->v.nr = 1;
723 sel->v.nalloc = -1;
724 return true;
727 /*! \brief
728 * Parses the values for a parameter that takes a constant number of values.
730 * \param[in] values List of values.
731 * \param param Parameter to parse.
732 * \param root Selection element to which child expressions are added.
733 * \param[in] scanner Scanner data structure.
734 * \returns true if the values were parsed successfully, false otherwise.
736 * For integer ranges, the sequence of numbers from the first to second value
737 * is stored, each as a separate value.
739 static bool
740 parse_values_std(const SelectionParserValueList &values,
741 gmx_ana_selparam_t *param,
742 const SelectionTreeElementPointer &root, void *scanner)
744 int i, j;
745 bool bDynamic;
747 /* Handle atom-valued parameters */
748 if (param->flags & SPAR_ATOMVAL)
750 if (values.size() > 1)
752 _gmx_selparser_error(scanner, "more than one value not supported");
753 return false;
755 if (values.front().hasExpressionValue())
757 SelectionTreeElementPointer child
758 = add_child(root, param, values.front().expr, scanner);
759 if (!child)
761 return false;
763 child->flags |= SEL_ALLOCVAL;
764 if (child->v.type != GROUP_VALUE && (child->flags & SEL_ATOMVAL))
766 /* Rest of the initialization is done during compilation in
767 * init_method(). */
768 /* TODO: Positions are not correctly handled */
769 param->val.nr = -1;
770 if (param->nvalptr)
772 *param->nvalptr = -1;
774 return true;
776 param->flags &= ~SPAR_ATOMVAL;
777 param->val.nr = 1;
778 if (param->nvalptr)
780 *param->nvalptr = 1;
782 param->nvalptr = NULL;
783 if (param->val.type == INT_VALUE || param->val.type == REAL_VALUE
784 || param->val.type == STR_VALUE)
786 _gmx_selvalue_reserve(&param->val, 1);
788 return set_expr_value_store(child, param, 0, scanner);
790 /* If we reach here, proceed with normal parameter handling */
791 param->val.nr = 1;
792 if (param->val.type == INT_VALUE || param->val.type == REAL_VALUE
793 || param->val.type == STR_VALUE)
795 _gmx_selvalue_reserve(&param->val, 1);
797 param->flags &= ~SPAR_ATOMVAL;
798 param->flags &= ~SPAR_DYNAMIC;
801 i = 0;
802 bDynamic = false;
803 SelectionParserValueList::const_iterator value;
804 for (value = values.begin(); value != values.end() && i < param->val.nr; ++value)
806 if (value->type != param->val.type)
808 _gmx_selparser_error(scanner, "incorrect value skipped");
809 continue;
811 if (value->hasExpressionValue())
813 SelectionTreeElementPointer child
814 = add_child(root, param, value->expr, scanner);
815 /* Check that the expression is valid */
816 if (!child)
818 return false;
820 if (!set_expr_value_store(child, param, i, scanner))
822 return false;
824 if (child->flags & SEL_DYNAMIC)
826 bDynamic = true;
829 else
831 /* Value is not an expression */
832 switch (value->type)
834 case INT_VALUE:
835 if (value->u.i.i1 <= value->u.i.i2)
837 for (j = value->u.i.i1; j <= value->u.i.i2 && i < param->val.nr; ++j)
839 param->val.u.i[i++] = j;
841 if (j != value->u.i.i2 + 1)
843 _gmx_selparser_error(scanner, "extra values skipped");
846 else
848 for (j = value->u.i.i1; j >= value->u.i.i2 && i < param->val.nr; --j)
850 param->val.u.i[i++] = j;
852 if (j != value->u.i.i2 - 1)
854 _gmx_selparser_error(scanner, "extra values skipped");
857 --i;
858 break;
859 case REAL_VALUE:
860 if (value->u.r.r1 != value->u.r.r2)
862 _gmx_selparser_error(scanner, "real ranges not supported");
863 return false;
865 param->val.u.r[i] = value->u.r.r1;
866 break;
867 case STR_VALUE:
868 param->val.u.s[i] = strdup(value->stringValue().c_str());
869 break;
870 case POS_VALUE:
871 gmx_ana_pos_init_const(&param->val.u.p[i], value->u.x);
872 break;
873 case NO_VALUE:
874 case GROUP_VALUE:
875 GMX_THROW(gmx::InternalError("Invalid non-expression value type"));
878 ++i;
880 if (value != values.end())
882 _gmx_selparser_error(scanner, "extra values");
883 return false;
885 if (i < param->val.nr)
887 _gmx_selparser_error(scanner, "not enough values");
888 return false;
890 if (!bDynamic)
892 param->flags &= ~SPAR_DYNAMIC;
894 if (param->nvalptr)
896 *param->nvalptr = param->val.nr;
898 param->nvalptr = NULL;
900 return true;
903 /*! \brief
904 * Parses the values for a boolean parameter.
906 * \param[in] name Name by which the parameter was given.
907 * \param[in] values List of values.
908 * \param param Parameter to parse.
909 * \param[in] scanner Scanner data structure.
910 * \returns true if the values were parsed successfully, false otherwise.
912 static bool
913 parse_values_bool(const std::string &name,
914 const SelectionParserValueList &values,
915 gmx_ana_selparam_t *param, void *scanner)
917 GMX_ASSERT(param->val.type == NO_VALUE,
918 "Boolean parser called for non-boolean parameter");
919 if (values.size() > 1 || (!values.empty() && values.front().type != INT_VALUE))
921 _gmx_selparser_error(scanner, "parameter takes only a yes/no/on/off/0/1 value");
922 return false;
925 bool bSetNo = false;
926 /* Check if the parameter name is given with a 'no' prefix */
927 if (name.length() > 2 && name[0] == 'n' && name[1] == 'o'
928 && name.compare(2, name.length() - 2, param->name) == 0)
930 bSetNo = true;
932 if (bSetNo && !values.empty())
934 _gmx_selparser_error(scanner, "parameter 'no%s' should not have a value",
935 param->name);
936 return false;
938 if (!values.empty() && values.front().u.i.i1 == 0)
940 bSetNo = true;
943 *param->val.u.b = bSetNo ? false : true;
944 return true;
947 /*! \brief
948 * Parses the values for an enumeration parameter.
950 * \param[in] values List of values.
951 * \param param Parameter to parse.
952 * \param[in] scanner Scanner data structure.
953 * \returns true if the values were parsed successfully, false otherwise.
955 static bool
956 parse_values_enum(const SelectionParserValueList &values,
957 gmx_ana_selparam_t *param,
958 void *scanner)
960 GMX_ASSERT(param->val.type == STR_VALUE,
961 "Enum parser called for non-string parameter");
962 if (values.size() != 1)
964 _gmx_selparser_error(scanner, "a single value is required");
965 return false;
967 const SelectionParserValue &value = values.front();
968 GMX_RELEASE_ASSERT(value.type == param->val.type,
969 "Invalid value type (should have been caught earlier)");
970 if (value.hasExpressionValue())
972 _gmx_selparser_error(scanner, "expression value for enumerated parameter not supported");
973 return false;
976 const std::string &svalue = value.stringValue();
977 int i = 1;
978 int match = 0;
979 while (param->val.u.s[i] != NULL)
981 if (gmx::startsWith(param->val.u.s[i], svalue))
983 /* Check if there is a duplicate match */
984 if (match > 0)
986 _gmx_selparser_error(scanner, "ambiguous value");
987 return false;
989 match = i;
991 ++i;
993 if (match == 0)
995 _gmx_selparser_error(scanner, "invalid value");
996 return false;
998 param->val.u.s[0] = param->val.u.s[match];
999 return true;
1002 /*! \brief
1003 * Replaces constant expressions with their values.
1005 * \param[in,out] values First element in the value list to process.
1007 static void
1008 convert_const_values(SelectionParserValueList *values)
1010 SelectionParserValueList::iterator value;
1011 for (value = values->begin(); value != values->end(); ++value)
1013 if (value->hasExpressionValue() && value->expr->v.type != GROUP_VALUE &&
1014 value->expr->type == SEL_CONST)
1016 SelectionTreeElementPointer expr = value->expr;
1017 switch (expr->v.type)
1019 case INT_VALUE:
1020 *value = SelectionParserValue::createInteger(expr->v.u.i[0]);
1021 break;
1022 case REAL_VALUE:
1023 *value = SelectionParserValue::createReal(expr->v.u.r[0]);
1024 break;
1025 case STR_VALUE:
1026 *value = SelectionParserValue::createString(expr->v.u.s[0]);
1027 break;
1028 case POS_VALUE:
1029 *value = SelectionParserValue::createPosition(expr->v.u.p->x[0]);
1030 break;
1031 default:
1032 GMX_THROW(gmx::InternalError(
1033 "Unsupported constant expression value type"));
1040 * \param pparams List of parameters from the selection parser.
1041 * \param[in] nparam Number of parameters in \p params.
1042 * \param params Array of parameters to parse.
1043 * \param root Selection element to which child expressions are added.
1044 * \param[in] scanner Scanner data structure.
1045 * \returns true if the parameters were parsed successfully, false otherwise.
1047 * Initializes the \p params array based on the parameters in \p pparams.
1048 * See the documentation of \c gmx_ana_selparam_t for different options
1049 * available for parsing.
1051 * The list \p pparams and any associated values are freed after the parameters
1052 * have been processed, no matter is there was an error or not.
1054 bool
1055 _gmx_sel_parse_params(const gmx::SelectionParserParameterList &pparams,
1056 int nparam, gmx_ana_selparam_t *params,
1057 const gmx::SelectionTreeElementPointer &root,
1058 void *scanner)
1060 gmx::MessageStringCollector *errors = _gmx_sel_lexer_error_reporter(scanner);
1061 gmx_ana_selparam_t *oparam;
1062 bool bOk, rc;
1063 int i;
1065 /* Check that the value pointers of SPAR_VARNUM parameters are NULL and
1066 * that they are not NULL for other parameters */
1067 bOk = true;
1068 for (i = 0; i < nparam; ++i)
1070 std::string contextStr = gmx::formatString("In parameter '%s'", params[i].name);
1071 gmx::MessageStringContext context(errors, contextStr);
1072 if (params[i].val.type != POS_VALUE && (params[i].flags & (SPAR_VARNUM | SPAR_ATOMVAL)))
1074 if (params[i].val.u.ptr != NULL)
1076 _gmx_selparser_error(scanner, "value pointer is not NULL "
1077 "although it should be for SPAR_VARNUM "
1078 "and SPAR_ATOMVAL parameters");
1080 if ((params[i].flags & SPAR_VARNUM)
1081 && (params[i].flags & SPAR_DYNAMIC) && !params[i].nvalptr)
1083 _gmx_selparser_error(scanner, "nvalptr is NULL but both "
1084 "SPAR_VARNUM and SPAR_DYNAMIC are specified");
1085 bOk = false;
1088 else
1090 if (params[i].val.u.ptr == NULL)
1092 _gmx_selparser_error(scanner, "value pointer is NULL");
1093 bOk = false;
1097 if (!bOk)
1099 return false;
1101 /* Parse the parameters */
1102 i = 0;
1103 SelectionParserParameterList::const_iterator pparam;
1104 for (pparam = pparams.begin(); pparam != pparams.end(); ++pparam)
1106 std::string contextStr;
1107 /* Find the parameter and make some checks */
1108 if (!pparam->name().empty())
1110 contextStr = gmx::formatString("In parameter '%s'", pparam->name().c_str());
1111 i = -1;
1112 oparam = gmx_ana_selparam_find(pparam->name().c_str(), nparam, params);
1114 else if (i >= 0)
1116 contextStr = gmx::formatString("In value %d", i + 1);
1117 oparam = &params[i];
1118 if (oparam->name != NULL)
1120 oparam = NULL;
1121 _gmx_selparser_error(scanner, "too many NULL parameters provided");
1122 bOk = false;
1123 continue;
1125 ++i;
1127 else
1129 _gmx_selparser_error(scanner, "all NULL parameters should appear in the beginning of the list");
1130 bOk = false;
1131 continue;
1133 gmx::MessageStringContext context(errors, contextStr);
1134 if (!oparam)
1136 _gmx_selparser_error(scanner, "unknown parameter skipped");
1137 bOk = false;
1138 continue;
1140 if (oparam->val.type != NO_VALUE && pparam->values().empty())
1142 _gmx_selparser_error(scanner, "no value provided");
1143 bOk = false;
1144 continue;
1146 if (oparam->flags & SPAR_SET)
1148 _gmx_selparser_error(scanner, "parameter set multiple times, extra values skipped");
1149 bOk = false;
1150 continue;
1152 oparam->flags |= SPAR_SET;
1153 /* Process the values for the parameter */
1154 convert_const_values(pparam->values_.get());
1155 if (convert_values(pparam->values_.get(), oparam->val.type, scanner) != 0)
1157 _gmx_selparser_error(scanner, "invalid value");
1158 bOk = false;
1159 continue;
1161 if (oparam->val.type == NO_VALUE)
1163 rc = parse_values_bool(pparam->name(), pparam->values(), oparam, scanner);
1165 else if (oparam->flags & SPAR_RANGES)
1167 rc = parse_values_range(pparam->values(), oparam, scanner);
1169 else if (oparam->flags & SPAR_VARNUM)
1171 if (pparam->values().size() == 1
1172 && pparam->values().front().hasExpressionValue())
1174 rc = parse_values_varnum_expr(pparam->values(), oparam, root, scanner);
1176 else
1178 rc = parse_values_varnum(pparam->values(), oparam, root, scanner);
1181 else if (oparam->flags & SPAR_ENUMVAL)
1183 rc = parse_values_enum(pparam->values(), oparam, scanner);
1185 else
1187 rc = parse_values_std(pparam->values(), oparam, root, scanner);
1189 if (!rc)
1191 bOk = false;
1194 /* Check that all required parameters are present */
1195 for (i = 0; i < nparam; ++i)
1197 if (!(params[i].flags & SPAR_OPTIONAL) && !(params[i].flags & SPAR_SET))
1199 _gmx_selparser_error(scanner, "required parameter '%s' not specified", params[i].name);
1200 bOk = false;
1204 return bOk;