C++ classes for selection parser values.
[gromacs.git] / src / gromacs / selection / params.cpp
blobf05c406ed9eca62c911a7a9fe5697bec1c60e0d5
1 /*
3 * This source code is part of
5 * G R O M A C S
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
31 /*! \internal \file
32 * \brief
33 * Implements functions in selparam.h.
35 * \author Teemu Murtola <teemu.murtola@cbr.su.se>
36 * \ingroup module_selection
38 #include <algorithm>
39 #include <string>
41 #include "gromacs/legacyheaders/smalloc.h"
42 #include "gromacs/legacyheaders/string2.h"
43 #include "gromacs/legacyheaders/vec.h"
45 #include "gromacs/selection/position.h"
46 #include "gromacs/selection/selmethod.h"
47 #include "gromacs/selection/selparam.h"
48 #include "gromacs/utility/errorcodes.h"
49 #include "gromacs/utility/exceptions.h"
50 #include "gromacs/utility/gmxassert.h"
51 #include "gromacs/utility/messagestringcollector.h"
52 #include "gromacs/utility/stringutil.h"
54 #include "parsetree.h"
55 #include "position.h"
56 #include "scanner.h"
57 #include "selelem.h"
59 using gmx::SelectionParserValue;
60 using gmx::SelectionParserValueList;
61 using gmx::SelectionParserParameterPointer;
62 using gmx::SelectionParserParameterList;
63 using gmx::SelectionTreeElement;
64 using gmx::SelectionTreeElementPointer;
66 /*!
67 * \param[in] name Name of the parameter to search.
68 * \param[in] nparam Number of parameters in the \p param array.
69 * \param[in] param Parameter array to search.
70 * \returns Pointer to the parameter in the \p param
71 * or NULL if no parameter with name \p name was found.
73 * The comparison is case-sensitive.
75 gmx_ana_selparam_t *
76 gmx_ana_selparam_find(const char *name, int nparam, gmx_ana_selparam_t *param)
78 int i;
80 if (nparam == 0)
82 return NULL;
84 /* Find the first non-null parameter */
85 i = 0;
86 while (i < nparam && param[i].name == NULL)
88 ++i;
90 /* Process the special case of a NULL parameter */
91 if (name == NULL)
93 return (i == 0) ? NULL : &param[i-1];
95 for ( ; i < nparam; ++i)
97 if (!strcmp(param[i].name, name))
99 return &param[i];
101 /* Check for 'no' prefix on boolean parameters */
102 if (param[i].val.type == NO_VALUE
103 && strlen(name) > 2 && name[0] == 'n' && name[1] == 'o'
104 && !strcmp(param[i].name, name+2))
106 return &param[i];
109 return NULL;
112 /*! \brief
113 * Does a type conversion on a SelectionParserValue.
115 * \param[in,out] value Value to convert.
116 * \param[in] type Type to convert to.
117 * \param[in] scanner Scanner data structure.
118 * \returns 0 on success, a non-zero value on error.
120 static int
121 convert_value(SelectionParserValue *value, e_selvalue_t type, void *scanner)
123 if (value->type == type || type == NO_VALUE)
125 return 0;
127 if (value->hasExpressionValue())
129 /* Conversion from atom selection to position using default
130 * reference positions. */
131 if (value->type == GROUP_VALUE && type == POS_VALUE)
133 SelectionTreeElementPointer expr =
134 _gmx_sel_init_position(value->expr, NULL, scanner);
135 // FIXME: Use exceptions
136 if (!expr)
138 return -1;
140 *value = SelectionParserValue::createExpr(expr);
141 return 0;
143 return -1;
145 else
147 /* Integers to floating point are easy */
148 if (value->type == INT_VALUE && type == REAL_VALUE)
150 *value = SelectionParserValue::createRealRange(value->u.i.i1,
151 value->u.i.i2);
152 return 0;
154 /* Reals that are integer-valued can also be converted */
155 if (value->type == REAL_VALUE && type == INT_VALUE)
157 int i1 = static_cast<int>(value->u.r.r1);
158 int i2 = static_cast<int>(value->u.r.r2);
159 if (gmx_within_tol(value->u.r.r1, i1, GMX_REAL_EPS)
160 && gmx_within_tol(value->u.r.r2, i2, GMX_REAL_EPS))
162 *value = SelectionParserValue::createIntegerRange(i1, i2);
163 return 0;
167 return -1;
170 /*! \brief
171 * Does a type conversion on a list of values.
173 * \param[in,out] values Values to convert.
174 * \param[in] type Type to convert to.
175 * \param[in] scanner Scanner data structure.
176 * \returns 0 on success, a non-zero value on error.
178 static int
179 convert_values(SelectionParserValueList *values, e_selvalue_t type, void *scanner)
181 int rc = 0;
182 SelectionParserValueList::iterator value;
183 for (value = values->begin(); value != values->end(); ++value)
185 int rc1 = convert_value(&*value, type, scanner);
186 if (rc1 != 0 && rc == 0)
188 rc = rc1;
191 /* FIXME: More informative error messages */
192 return rc;
195 /*! \brief
196 * Adds a child element for a parameter, keeping the parameter order.
198 * \param[in,out] root Root element to which the child is added.
199 * \param[in] child Child to add.
200 * \param[in] param Parameter for which this child is a value.
202 * Puts \p child in the child list of \p root such that the list remains
203 * in the same order as the corresponding parameters.
205 static void
206 place_child(const SelectionTreeElementPointer &root,
207 const SelectionTreeElementPointer &child,
208 gmx_ana_selparam_t *param)
210 gmx_ana_selparam_t *ps;
211 int n;
213 ps = root->u.expr.method->param;
214 n = param - ps;
215 /* Put the child element in the correct place */
216 if (!root->child || n < root->child->u.param - ps)
218 child->next = root->child;
219 root->child = child;
221 else
223 SelectionTreeElementPointer prev = root->child;
224 while (prev->next && prev->next->u.param - ps >= n)
226 prev = prev->next;
228 child->next = prev->next;
229 prev->next = child;
233 /*! \brief
234 * Comparison function for sorting integer ranges.
236 * \param[in] a Pointer to the first range.
237 * \param[in] b Pointer to the second range.
238 * \returns -1, 0, or 1 depending on the relative order of \p a and \p b.
240 * The ranges are primarily sorted based on their starting point, and
241 * secondarily based on length (longer ranges come first).
243 static int
244 cmp_int_range(const void *a, const void *b)
246 if (((int *)a)[0] < ((int *)b)[0])
248 return -1;
250 if (((int *)a)[0] > ((int *)b)[0])
252 return 1;
254 if (((int *)a)[1] > ((int *)b)[1])
256 return -1;
258 return 0;
261 /*! \brief
262 * Comparison function for sorting real ranges.
264 * \param[in] a Pointer to the first range.
265 * \param[in] b Pointer to the second range.
266 * \returns -1, 0, or 1 depending on the relative order of \p a and \p b.
268 * The ranges are primarily sorted based on their starting point, and
269 * secondarily based on length (longer ranges come first).
271 static int
272 cmp_real_range(const void *a, const void *b)
274 if (((real *)a)[0] < ((real *)b)[0])
276 return -1;
278 if (((real *)a)[0] > ((real *)b)[0])
280 return 1;
282 if (((real *)a)[1] > ((real *)b)[1])
284 return -1;
286 return 0;
289 /*! \brief
290 * Parses the values for a parameter that takes integer or real ranges.
292 * \param[in] values List of values.
293 * \param param Parameter to parse.
294 * \param[in] scanner Scanner data structure.
295 * \returns true if the values were parsed successfully, false otherwise.
297 static bool
298 parse_values_range(const SelectionParserValueList &values,
299 gmx_ana_selparam_t *param, void *scanner)
301 int *idata;
302 real *rdata;
303 int i, j, n;
305 param->flags &= ~SPAR_DYNAMIC;
306 if (param->val.type != INT_VALUE && param->val.type != REAL_VALUE)
308 GMX_ERROR_NORET(gmx::eeInternalError, "Invalid range parameter type");
309 return false;
311 idata = NULL;
312 rdata = NULL;
313 if (param->val.type == INT_VALUE)
315 snew(idata, values.size()*2);
317 else
319 snew(rdata, values.size()*2);
321 i = 0;
322 SelectionParserValueList::const_iterator value;
323 for (value = values.begin(); value != values.end(); ++value)
325 if (value->hasExpressionValue())
327 _gmx_selparser_error(scanner, "expressions not supported within range parameters");
328 return false;
330 if (value->type != param->val.type)
332 GMX_ERROR_NORET(gmx::eeInternalError, "Invalid range value type");
333 return false;
335 if (param->val.type == INT_VALUE)
337 int i1 = std::min(value->u.i.i1, value->u.i.i2);
338 int i2 = std::max(value->u.i.i1, value->u.i.i2);
339 /* Check if the new range overlaps or extends the previous one */
340 if (i > 0 && i1 <= idata[i-1]+1 && i2 >= idata[i-2]-1)
342 idata[i-2] = std::min(idata[i-2], i1);
343 idata[i-1] = std::max(idata[i-1], i2);
345 else
347 idata[i++] = i1;
348 idata[i++] = i2;
351 else
353 real r1 = std::min(value->u.r.r1, value->u.r.r2);
354 real r2 = std::max(value->u.r.r1, value->u.r.r2);
355 /* Check if the new range overlaps or extends the previous one */
356 if (i > 0 && r1 <= rdata[i-1] && r2 >= rdata[i-2])
358 rdata[i-2] = std::min(rdata[i-2], r1);
359 rdata[i-1] = std::max(rdata[i-1], r2);
361 else
363 rdata[i++] = r1;
364 rdata[i++] = r2;
368 n = i/2;
369 /* Sort the ranges and merge consequent ones */
370 if (param->val.type == INT_VALUE)
372 qsort(idata, n, 2*sizeof(int), &cmp_int_range);
373 for (i = j = 2; i < 2*n; i += 2)
375 if (idata[j-1]+1 >= idata[i])
377 if (idata[i+1] > idata[j-1])
379 idata[j-1] = idata[i+1];
382 else
384 idata[j] = idata[i];
385 idata[j+1] = idata[i+1];
386 j += 2;
390 else
392 qsort(rdata, n, 2*sizeof(real), &cmp_real_range);
393 for (i = j = 2; i < 2*n; i += 2)
395 if (rdata[j-1]+1 >= rdata[i])
397 if (rdata[i+1] > rdata[j-1])
399 rdata[j-1] = rdata[i+1];
402 else
404 rdata[j] = rdata[i];
405 rdata[j+1] = rdata[i+1];
406 j += 2;
410 n = j/2;
411 /* Store the values */
412 if (param->flags & SPAR_VARNUM)
414 param->val.nr = n;
415 if (param->val.type == INT_VALUE)
417 srenew(idata, j);
418 _gmx_selvalue_setstore_alloc(&param->val, idata, j);
420 else
422 srenew(rdata, j);
423 _gmx_selvalue_setstore_alloc(&param->val, rdata, j);
426 else
428 if (n != param->val.nr)
430 _gmx_selparser_error(scanner, "the value should consist of exactly one range");
431 sfree(idata);
432 sfree(rdata);
433 return false;
435 if (param->val.type == INT_VALUE)
437 memcpy(param->val.u.i, idata, 2*n*sizeof(int));
438 sfree(idata);
440 else
442 memcpy(param->val.u.r, rdata, 2*n*sizeof(real));
443 sfree(rdata);
446 if (param->nvalptr)
448 *param->nvalptr = param->val.nr;
450 param->nvalptr = NULL;
452 return true;
455 /*! \brief
456 * Parses the values for a parameter that takes a variable number of values.
458 * \param[in] values List of values.
459 * \param param Parameter to parse.
460 * \param root Selection element to which child expressions are added.
461 * \param[in] scanner Scanner data structure.
462 * \returns true if the values were parsed successfully, false otherwise.
464 * For integer ranges, the sequence of numbers from the first to second value
465 * is stored, each as a separate value.
467 static bool
468 parse_values_varnum(const SelectionParserValueList &values,
469 gmx_ana_selparam_t *param,
470 const SelectionTreeElementPointer &root,
471 void *scanner)
473 int i, j;
475 param->flags &= ~SPAR_DYNAMIC;
476 /* Compute number of values, considering also integer ranges. */
477 size_t valueCount = values.size();
478 if (param->val.type == INT_VALUE)
480 SelectionParserValueList::const_iterator value;
481 for (value = values.begin(); value != values.end(); ++value)
483 if (value->type == INT_VALUE && !value->hasExpressionValue())
485 valueCount += abs(value->u.i.i2 - value->u.i.i1);
490 /* Check that the value type is actually implemented */
491 if (param->val.type != INT_VALUE && param->val.type != REAL_VALUE
492 && param->val.type != STR_VALUE && param->val.type != POS_VALUE)
494 GMX_ERROR_NORET(gmx::eeInternalError,
495 "Variable-count value type not implemented");
496 return false;
499 /* Reserve appropriate amount of memory */
500 if (param->val.type == POS_VALUE)
502 gmx_ana_pos_reserve(param->val.u.p, valueCount, 0);
503 gmx_ana_pos_set_nr(param->val.u.p, valueCount);
504 gmx_ana_indexmap_init(&param->val.u.p->m, NULL, NULL, INDEX_UNKNOWN);
506 else
508 _gmx_selvalue_reserve(&param->val, valueCount);
511 i = 0;
512 SelectionParserValueList::const_iterator value;
513 for (value = values.begin(); value != values.end(); ++value)
515 if (value->hasExpressionValue())
517 _gmx_selparser_error(scanner, "expressions not supported within value lists");
518 return false;
520 if (value->type != param->val.type)
522 GMX_ERROR_NORET(gmx::eeInternalError, "Invalid value type");
523 return false;
525 switch (param->val.type)
527 case INT_VALUE:
528 if (value->u.i.i1 <= value->u.i.i2)
530 for (j = value->u.i.i1; j <= value->u.i.i2; ++j)
532 param->val.u.i[i++] = j;
535 else
537 for (j = value->u.i.i1; j >= value->u.i.i2; --j)
539 param->val.u.i[i++] = j;
542 break;
543 case REAL_VALUE:
544 if (value->u.r.r1 != value->u.r.r2)
546 _gmx_selparser_error(scanner, "real ranges not supported");
547 return false;
549 param->val.u.r[i++] = value->u.r.r1;
550 break;
551 case STR_VALUE:
552 param->val.u.s[i++] = strdup(value->stringValue().c_str());
553 break;
554 case POS_VALUE: copy_rvec(value->u.x, param->val.u.p->x[i++]); break;
555 default: /* Should not be reached */
556 GMX_ERROR_NORET(gmx::eeInternalError, "Invalid value type");
557 return false;
560 param->val.nr = i;
561 if (param->nvalptr)
563 *param->nvalptr = param->val.nr;
565 param->nvalptr = NULL;
566 /* Create a dummy child element to store the string values.
567 * This element is responsible for freeing the values, but carries no
568 * other function. */
569 if (param->val.type == STR_VALUE)
571 SelectionTreeElementPointer child(new SelectionTreeElement(SEL_CONST));
572 _gmx_selelem_set_vtype(child, STR_VALUE);
573 child->setName(param->name);
574 child->flags &= ~SEL_ALLOCVAL;
575 child->flags |= SEL_FLAGSSET | SEL_VARNUMVAL | SEL_ALLOCDATA;
576 child->v.nr = param->val.nr;
577 _gmx_selvalue_setstore(&child->v, param->val.u.s);
578 /* Because the child is not group-valued, the u union is not used
579 * for anything, so we can abuse it by storing the parameter value
580 * as place_child() expects, but this is really ugly... */
581 child->u.param = param;
582 place_child(root, child, param);
585 return true;
588 /*! \brief
589 * Adds a new subexpression reference to a selection element.
591 * \param[in,out] root Root element to which the subexpression is added.
592 * \param[in] param Parameter for which this expression is a value.
593 * \param[in] expr Expression to add.
594 * \param[in] scanner Scanner data structure.
595 * \returns The created child element.
597 * Creates a new \ref SEL_SUBEXPRREF element and adds it into the child
598 * list of \p root.
599 * If \p expr is already a \ref SEL_SUBEXPRREF, it is used as it is.
600 * \ref SEL_ALLOCVAL is cleared for the returned element.
602 static SelectionTreeElementPointer
603 add_child(const SelectionTreeElementPointer &root, gmx_ana_selparam_t *param,
604 const SelectionTreeElementPointer &expr, void *scanner)
606 GMX_RELEASE_ASSERT(root->type == SEL_EXPRESSION || root->type == SEL_MODIFIER,
607 "Unsupported root element for selection parameter parser");
608 SelectionTreeElementPointer child;
609 /* Create a subexpression reference element if necessary */
610 if (expr->type == SEL_SUBEXPRREF)
612 child = expr;
614 else
616 child.reset(new SelectionTreeElement(SEL_SUBEXPRREF));
617 _gmx_selelem_set_vtype(child, expr->v.type);
618 child->child = expr;
620 /* Setup the child element */
621 child->flags &= ~SEL_ALLOCVAL;
622 child->u.param = param;
623 if (child->v.type != param->val.type)
625 _gmx_selparser_error(scanner, "invalid expression value");
626 // FIXME: Use exceptions.
627 return SelectionTreeElementPointer();
629 _gmx_selelem_update_flags(child, scanner);
630 if ((child->flags & SEL_DYNAMIC) && !(param->flags & SPAR_DYNAMIC))
632 _gmx_selparser_error(scanner, "dynamic values not supported");
633 // FIXME: Use exceptions.
634 return SelectionTreeElementPointer();
636 if (!(child->flags & SEL_DYNAMIC))
638 param->flags &= ~SPAR_DYNAMIC;
640 /* Put the child element in the correct place */
641 place_child(root, child, param);
642 return child;
645 /*! \brief
646 * Parses an expression value for a parameter that takes a variable number of values.
648 * \param[in] values List of values.
649 * \param param Parameter to parse.
650 * \param root Selection element to which child expressions are added.
651 * \param[in] scanner Scanner data structure.
652 * \returns true if the values were parsed successfully, false otherwise.
654 static bool
655 parse_values_varnum_expr(const SelectionParserValueList &values,
656 gmx_ana_selparam_t *param,
657 const SelectionTreeElementPointer &root,
658 void *scanner)
660 if (values.size() != 1 || !values.front().hasExpressionValue())
662 GMX_ERROR_NORET(gmx::eeInternalError, "Invalid expression value");
663 return false;
666 SelectionTreeElementPointer child
667 = add_child(root, param, values.front().expr, scanner);
668 if (!child)
670 return false;
673 /* Process single-valued expressions */
674 /* TODO: We should also handle SEL_SINGLEVAL expressions here */
675 if (child->v.type == POS_VALUE || child->v.type == GROUP_VALUE)
677 /* Set the value storage */
678 _gmx_selvalue_setstore(&child->v, param->val.u.ptr);
679 param->val.nr = 1;
680 if (param->nvalptr)
682 *param->nvalptr = param->val.nr;
684 param->nvalptr = NULL;
685 return true;
688 if (!(child->flags & SEL_VARNUMVAL))
690 _gmx_selparser_error(scanner, "invalid expression value");
691 return false;
694 child->flags |= SEL_ALLOCVAL;
695 param->val.nr = -1;
696 *param->nvalptr = param->val.nr;
697 /* Rest of the initialization is done during compilation in
698 * init_method(). */
700 return true;
703 /*! \brief
704 * Initializes the storage of an expression value.
706 * \param[in,out] sel Selection element that evaluates the value.
707 * \param[in] param Parameter to receive the value.
708 * \param[in] i The value of \p sel evaluates the value \p i for
709 * \p param.
710 * \param[in] scanner Scanner data structure.
712 * Initializes the data pointer of \p sel such that the result is stored
713 * as the value \p i of \p param.
714 * This function is used internally by parse_values_std().
716 static bool
717 set_expr_value_store(const SelectionTreeElementPointer &sel,
718 gmx_ana_selparam_t *param, int i, void *scanner)
720 if (sel->v.type != GROUP_VALUE && !(sel->flags & SEL_SINGLEVAL))
722 _gmx_selparser_error(scanner, "invalid expression value");
723 return false;
725 switch (sel->v.type)
727 case INT_VALUE: sel->v.u.i = &param->val.u.i[i]; break;
728 case REAL_VALUE: sel->v.u.r = &param->val.u.r[i]; break;
729 case STR_VALUE: sel->v.u.s = &param->val.u.s[i]; break;
730 case POS_VALUE: sel->v.u.p = &param->val.u.p[i]; break;
731 case GROUP_VALUE: sel->v.u.g = &param->val.u.g[i]; break;
732 default: /* Error */
733 GMX_ERROR_NORET(gmx::eeInternalError, "Invalid value type");
734 return false;
736 sel->v.nr = 1;
737 sel->v.nalloc = -1;
738 return true;
741 /*! \brief
742 * Parses the values for a parameter that takes a constant number of values.
744 * \param[in] values List of values.
745 * \param param Parameter to parse.
746 * \param root Selection element to which child expressions are added.
747 * \param[in] scanner Scanner data structure.
748 * \returns true if the values were parsed successfully, false otherwise.
750 * For integer ranges, the sequence of numbers from the first to second value
751 * is stored, each as a separate value.
753 static bool
754 parse_values_std(const SelectionParserValueList &values,
755 gmx_ana_selparam_t *param,
756 const SelectionTreeElementPointer &root, void *scanner)
758 int i, j;
759 bool bDynamic;
761 /* Handle atom-valued parameters */
762 if (param->flags & SPAR_ATOMVAL)
764 if (values.size() > 1)
766 _gmx_selparser_error(scanner, "more than one value not supported");
767 return false;
769 if (values.front().hasExpressionValue())
771 SelectionTreeElementPointer child
772 = add_child(root, param, values.front().expr, scanner);
773 if (!child)
775 return false;
777 child->flags |= SEL_ALLOCVAL;
778 if (child->v.type != GROUP_VALUE && (child->flags & SEL_ATOMVAL))
780 /* Rest of the initialization is done during compilation in
781 * init_method(). */
782 /* TODO: Positions are not correctly handled */
783 param->val.nr = -1;
784 if (param->nvalptr)
786 *param->nvalptr = -1;
788 return true;
790 param->flags &= ~SPAR_ATOMVAL;
791 param->val.nr = 1;
792 if (param->nvalptr)
794 *param->nvalptr = 1;
796 param->nvalptr = NULL;
797 if (param->val.type == INT_VALUE || param->val.type == REAL_VALUE
798 || param->val.type == STR_VALUE)
800 _gmx_selvalue_reserve(&param->val, 1);
802 return set_expr_value_store(child, param, 0, scanner);
804 /* If we reach here, proceed with normal parameter handling */
805 param->val.nr = 1;
806 if (param->val.type == INT_VALUE || param->val.type == REAL_VALUE
807 || param->val.type == STR_VALUE)
809 _gmx_selvalue_reserve(&param->val, 1);
811 param->flags &= ~SPAR_ATOMVAL;
812 param->flags &= ~SPAR_DYNAMIC;
815 i = 0;
816 bDynamic = false;
817 SelectionParserValueList::const_iterator value;
818 for (value = values.begin(); value != values.end() && i < param->val.nr; ++value)
820 if (value->type != param->val.type)
822 _gmx_selparser_error(scanner, "incorrect value skipped");
823 continue;
825 if (value->hasExpressionValue())
827 SelectionTreeElementPointer child
828 = add_child(root, param, value->expr, scanner);
829 /* Check that the expression is valid */
830 if (!child)
832 return false;
834 if (!set_expr_value_store(child, param, i, scanner))
836 return false;
838 if (child->flags & SEL_DYNAMIC)
840 bDynamic = true;
843 else
845 /* Value is not an expression */
846 switch (value->type)
848 case INT_VALUE:
849 if (value->u.i.i1 <= value->u.i.i2)
851 for (j = value->u.i.i1; j <= value->u.i.i2 && i < param->val.nr; ++j)
853 param->val.u.i[i++] = j;
855 if (j != value->u.i.i2 + 1)
857 _gmx_selparser_error(scanner, "extra values skipped");
860 else
862 for (j = value->u.i.i1; j >= value->u.i.i2 && i < param->val.nr; --j)
864 param->val.u.i[i++] = j;
866 if (j != value->u.i.i2 - 1)
868 _gmx_selparser_error(scanner, "extra values skipped");
871 --i;
872 break;
873 case REAL_VALUE:
874 if (value->u.r.r1 != value->u.r.r2)
876 _gmx_selparser_error(scanner, "real ranges not supported");
877 return false;
879 param->val.u.r[i] = value->u.r.r1;
880 break;
881 case STR_VALUE:
882 param->val.u.s[i] = strdup(value->stringValue().c_str());
883 break;
884 case POS_VALUE:
885 gmx_ana_pos_init_const(&param->val.u.p[i], value->u.x);
886 break;
887 case NO_VALUE:
888 case GROUP_VALUE:
889 GMX_ERROR_NORET(gmx::eeInternalError,
890 "Invalid non-expression value");
891 return false;
894 ++i;
896 if (value != values.end())
898 _gmx_selparser_error(scanner, "extra values");
899 return false;
901 if (i < param->val.nr)
903 _gmx_selparser_error(scanner, "not enough values");
904 return false;
906 if (!bDynamic)
908 param->flags &= ~SPAR_DYNAMIC;
910 if (param->nvalptr)
912 *param->nvalptr = param->val.nr;
914 param->nvalptr = NULL;
916 return true;
919 /*! \brief
920 * Parses the values for a boolean parameter.
922 * \param[in] name Name by which the parameter was given.
923 * \param[in] values List of values.
924 * \param param Parameter to parse.
925 * \param[in] scanner Scanner data structure.
926 * \returns true if the values were parsed successfully, false otherwise.
928 static bool
929 parse_values_bool(const std::string &name,
930 const SelectionParserValueList &values,
931 gmx_ana_selparam_t *param, void *scanner)
933 bool bSetNo;
935 if (param->val.type != NO_VALUE)
937 GMX_ERROR_NORET(gmx::eeInternalError, "Invalid boolean parameter");
938 return false;
940 if (values.size() > 1 || (!values.empty() && values.front().type != INT_VALUE))
942 _gmx_selparser_error(scanner, "parameter takes only a yes/no/on/off/0/1 value");
943 return false;
946 bSetNo = false;
947 /* Check if the parameter name is given with a 'no' prefix */
948 if (name.length() > 2 && name[0] == 'n' && name[1] == 'o'
949 && name.compare(2, name.length() - 2, param->name) == 0)
951 bSetNo = true;
953 if (bSetNo && !values.empty())
955 _gmx_selparser_error(scanner, "parameter 'no%s' should not have a value",
956 param->name);
957 return false;
959 if (!values.empty() && values.front().u.i.i1 == 0)
961 bSetNo = true;
964 *param->val.u.b = bSetNo ? false : true;
965 return true;
968 /*! \brief
969 * Parses the values for an enumeration parameter.
971 * \param[in] values List of values.
972 * \param param Parameter to parse.
973 * \param[in] scanner Scanner data structure.
974 * \returns true if the values were parsed successfully, false otherwise.
976 static bool
977 parse_values_enum(const SelectionParserValueList &values,
978 gmx_ana_selparam_t *param,
979 void *scanner)
981 if (values.size() != 1)
983 _gmx_selparser_error(scanner, "a single value is required");
984 return false;
986 const SelectionParserValue &value = values.front();
987 if (value.type != STR_VALUE || param->val.type != STR_VALUE)
989 GMX_ERROR_NORET(gmx::eeInternalError, "Invalid enum parameter");
990 return false;
992 if (value.hasExpressionValue())
994 _gmx_selparser_error(scanner, "expression value for enumerated parameter not supported");
995 return false;
998 const std::string &svalue = value.stringValue();
999 int i = 1;
1000 int match = 0;
1001 while (param->val.u.s[i] != NULL)
1003 if (gmx::startsWith(param->val.u.s[i], svalue))
1005 /* Check if there is a duplicate match */
1006 if (match > 0)
1008 _gmx_selparser_error(scanner, "ambiguous value");
1009 return false;
1011 match = i;
1013 ++i;
1015 if (match == 0)
1017 _gmx_selparser_error(scanner, "invalid value");
1018 return false;
1020 param->val.u.s[0] = param->val.u.s[match];
1021 return true;
1024 /*! \brief
1025 * Replaces constant expressions with their values.
1027 * \param[in,out] values First element in the value list to process.
1029 static void
1030 convert_const_values(SelectionParserValueList *values)
1032 SelectionParserValueList::iterator value;
1033 for (value = values->begin(); value != values->end(); ++value)
1035 if (value->hasExpressionValue() && value->expr->v.type != GROUP_VALUE &&
1036 value->expr->type == SEL_CONST)
1038 SelectionTreeElementPointer expr = value->expr;
1039 switch (expr->v.type)
1041 case INT_VALUE:
1042 *value = SelectionParserValue::createInteger(expr->v.u.i[0]);
1043 break;
1044 case REAL_VALUE:
1045 *value = SelectionParserValue::createReal(expr->v.u.r[0]);
1046 break;
1047 case STR_VALUE:
1048 *value = SelectionParserValue::createString(expr->v.u.s[0]);
1049 break;
1050 case POS_VALUE:
1051 *value = SelectionParserValue::createPosition(expr->v.u.p->x[0]);
1052 break;
1053 default:
1054 GMX_ERROR_NORET(gmx::eeInternalError,
1055 "Unsupported value type");
1056 break;
1063 * \param pparams List of parameters from the selection parser.
1064 * \param[in] nparam Number of parameters in \p params.
1065 * \param params Array of parameters to parse.
1066 * \param root Selection element to which child expressions are added.
1067 * \param[in] scanner Scanner data structure.
1068 * \returns true if the parameters were parsed successfully, false otherwise.
1070 * Initializes the \p params array based on the parameters in \p pparams.
1071 * See the documentation of \c gmx_ana_selparam_t for different options
1072 * available for parsing.
1074 * The list \p pparams and any associated values are freed after the parameters
1075 * have been processed, no matter is there was an error or not.
1077 bool
1078 _gmx_sel_parse_params(const SelectionParserParameterList &pparams,
1079 int nparam, gmx_ana_selparam_t *params,
1080 const SelectionTreeElementPointer &root, void *scanner)
1082 gmx::MessageStringCollector *errors = _gmx_sel_lexer_error_reporter(scanner);
1083 gmx_ana_selparam_t *oparam;
1084 bool bOk, rc;
1085 int i;
1087 /* Check that the value pointers of SPAR_VARNUM parameters are NULL and
1088 * that they are not NULL for other parameters */
1089 bOk = true;
1090 for (i = 0; i < nparam; ++i)
1092 std::string contextStr = gmx::formatString("In parameter '%s'", params[i].name);
1093 gmx::MessageStringContext context(errors, contextStr);
1094 if (params[i].val.type != POS_VALUE && (params[i].flags & (SPAR_VARNUM | SPAR_ATOMVAL)))
1096 if (params[i].val.u.ptr != NULL)
1098 _gmx_selparser_error(scanner, "value pointer is not NULL "
1099 "although it should be for SPAR_VARNUM "
1100 "and SPAR_ATOMVAL parameters");
1102 if ((params[i].flags & SPAR_VARNUM)
1103 && (params[i].flags & SPAR_DYNAMIC) && !params[i].nvalptr)
1105 _gmx_selparser_error(scanner, "nvalptr is NULL but both "
1106 "SPAR_VARNUM and SPAR_DYNAMIC are specified");
1107 bOk = false;
1110 else
1112 if (params[i].val.u.ptr == NULL)
1114 _gmx_selparser_error(scanner, "value pointer is NULL");
1115 bOk = false;
1119 if (!bOk)
1121 return false;
1123 /* Parse the parameters */
1124 i = 0;
1125 SelectionParserParameterList::const_iterator piter;
1126 for (piter = pparams.begin(); piter != pparams.end(); ++piter)
1128 const SelectionParserParameterPointer &pparam = *piter;
1129 std::string contextStr;
1130 /* Find the parameter and make some checks */
1131 if (!pparam->name().empty())
1133 contextStr = gmx::formatString("In parameter '%s'", pparam->name().c_str());
1134 i = -1;
1135 oparam = gmx_ana_selparam_find(pparam->name().c_str(), nparam, params);
1137 else if (i >= 0)
1139 contextStr = gmx::formatString("In value %d", i + 1);
1140 oparam = &params[i];
1141 if (oparam->name != NULL)
1143 oparam = NULL;
1144 _gmx_selparser_error(scanner, "too many NULL parameters provided");
1145 bOk = false;
1146 continue;
1148 ++i;
1150 else
1152 _gmx_selparser_error(scanner, "all NULL parameters should appear in the beginning of the list");
1153 bOk = false;
1154 continue;
1156 gmx::MessageStringContext context(errors, contextStr);
1157 if (!oparam)
1159 _gmx_selparser_error(scanner, "unknown parameter skipped");
1160 bOk = false;
1161 continue;
1163 if (oparam->val.type != NO_VALUE && pparam->values().empty())
1165 _gmx_selparser_error(scanner, "no value provided");
1166 bOk = false;
1167 continue;
1169 if (oparam->flags & SPAR_SET)
1171 _gmx_selparser_error(scanner, "parameter set multiple times, extra values skipped");
1172 bOk = false;
1173 continue;
1175 oparam->flags |= SPAR_SET;
1176 /* Process the values for the parameter */
1177 convert_const_values(pparam->values_.get());
1178 if (convert_values(pparam->values_.get(), oparam->val.type, scanner) != 0)
1180 _gmx_selparser_error(scanner, "invalid value");
1181 bOk = false;
1182 continue;
1184 if (oparam->val.type == NO_VALUE)
1186 rc = parse_values_bool(pparam->name(), pparam->values(), oparam, scanner);
1188 else if (oparam->flags & SPAR_RANGES)
1190 rc = parse_values_range(pparam->values(), oparam, scanner);
1192 else if (oparam->flags & SPAR_VARNUM)
1194 if (pparam->values().size() == 1
1195 && pparam->values().front().hasExpressionValue())
1197 rc = parse_values_varnum_expr(pparam->values(), oparam, root, scanner);
1199 else
1201 rc = parse_values_varnum(pparam->values(), oparam, root, scanner);
1204 else if (oparam->flags & SPAR_ENUMVAL)
1206 rc = parse_values_enum(pparam->values(), oparam, scanner);
1208 else
1210 rc = parse_values_std(pparam->values(), oparam, root, scanner);
1212 if (!rc)
1214 bOk = false;
1217 /* Check that all required parameters are present */
1218 for (i = 0; i < nparam; ++i)
1220 if (!(params[i].flags & SPAR_OPTIONAL) && !(params[i].flags & SPAR_SET))
1222 _gmx_selparser_error(scanner, "required parameter '%s' not specified", params[i].name);
1223 bOk = false;
1227 return bOk;