Remove some hardcoded selection string buffers
[gromacs.git] / src / gromacs / selection / parsetree.cpp
blob3035805bd4e11ce0f19f0631faff4351c6952f91
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 parsetree.h.
39 * \author Teemu Murtola <teemu.murtola@gmail.com>
40 * \ingroup module_selection
42 /*! \internal
43 * \page page_module_selection_parser Selection parsing
45 * The selection parser is implemented in the following files:
46 * - scanner.l:
47 * Tokenizer implemented using Flex, splits the input into tokens
48 * (scanner.c and scanner_flex.h are generated from this file).
49 * - scanner.h, scanner_internal.h, scanner_internal.cpp:
50 * Helper functions for scanner.l and for interfacing between
51 * scanner.l and parser.y. Functions in scanner_internal.h are only
52 * used from scanner.l, while scanner.h is used from the parser.
53 * - symrec.h, symrec.cpp:
54 * Functions used by the tokenizer to handle the symbol table, i.e.,
55 * the recognized keywords. Some basic keywords are hardcoded into
56 * scanner.l, but all method and variable references go through the
57 * symbol table, as do position evaluation keywords.
58 * - parser.y:
59 * Semantic rules for parsing the grammar
60 * (parser.cpp and parser.h are generated from this file by Bison).
61 * - parsetree.h, parsetree.cpp:
62 * Functions called from actions in parser.y to construct the
63 * evaluation elements corresponding to different grammar elements.
64 * - params.c:
65 * Defines a function that processes the parameters of selection
66 * methods and initializes the children of the method element.
67 * - selectioncollection.h, selectioncollection.cpp:
68 * These files define the high-level public interface to the parser
69 * through SelectionCollection::parseFromStdin(),
70 * SelectionCollection::parseFromFile() and
71 * SelectionCollection::parseFromString().
73 * The basic control flow in the parser is as follows: when a parser function
74 * in SelectionCollection gets called, it performs some
75 * initialization, and then calls the _gmx_sel_yyparse() function generated
76 * by Bison. This function then calls _gmx_sel_yylex() to repeatedly read
77 * tokens from the input (more complex tasks related to token recognition
78 * and bookkeeping are done by functions in scanner_internal.cpp) and uses the
79 * grammar rules to decide what to do with them. Whenever a grammar rule
80 * matches, a corresponding function in parsetree.cpp is called to construct
81 * either a temporary representation for the object or a
82 * gmx::SelectionTreeElement object
83 * (some simple rules are handled internally in parser.y).
84 * When a complete selection has been parsed, the functions in parsetree.cpp
85 * also take care of updating the ::gmx_ana_selcollection_t structure
86 * appropriately.
88 * The rest of this page describes the resulting gmx::SelectionTreeElement
89 * object tree.
90 * Before the selections can be evaluated, this tree needs to be passed to
91 * the selection compiler, which is described on a separate page:
92 * \ref page_module_selection_compiler
95 * \section selparser_tree Element tree constructed by the parser
97 * The parser initializes the following fields in all selection elements:
98 * gmx::SelectionTreeElement::name, gmx::SelectionTreeElement::type,
99 * gmx::SelectionTreeElement::v\c .type,
100 * gmx::SelectionTreeElement::flags, gmx::SelectionTreeElement::child, and
101 * gmx::SelectionTreeElement::next.
102 * Some other fields are also initialized for particular element types as
103 * discussed below.
104 * Fields that are not initialized are set to zero, NULL, or other similar
105 * value.
108 * \subsection selparser_tree_root Root elements
110 * The parser creates a \ref SEL_ROOT selection element for each variable
111 * assignment and each selection. However, there are two exceptions that do
112 * not result in a \ref SEL_ROOT element (in these cases, only the symbol
113 * table is modified):
114 * - Variable assignments that assign a variable to another variable.
115 * - Variable assignments that assign a non-group constant.
117 * The \ref SEL_ROOT elements are linked together in a chain in the same order
118 * as in the input.
120 * The children of the \ref SEL_ROOT elements can be used to distinguish
121 * the two types of root elements from each other:
122 * - For variable assignments, the first and only child is always
123 * a \ref SEL_SUBEXPR element.
124 * - For selections, the first child is a \ref SEL_EXPRESSION or a
125 * \ref SEL_MODIFIER element that evaluates the final positions (if the
126 * selection defines a constant position, the child is a \ref SEL_CONST).
127 * The rest of the children are \ref SEL_MODIFIER elements with
128 * \ref NO_VALUE, in the order given by the user.
130 * The name of the selection/variable is stored in
131 * gmx::SelectionTreeElement::cgrp\c .name.
132 * It is set to either the name provided by the user or the selection string
133 * for selections not explicitly named by the user.
134 * \ref SEL_ROOT or \ref SEL_SUBEXPR elements do not appear anywhere else.
137 * \subsection selparser_tree_const Constant elements
139 * \ref SEL_CONST elements are created for every constant that is required
140 * for later evaluation.
141 * Currently, \ref SEL_CONST elements can be present for
142 * - selections that consist of a constant position,
143 * - \ref GROUP_VALUE method parameters if provided using external index
144 * groups,
146 * For group-valued elements, the value is stored in
147 * gmx::SelectionTreeElement::cgrp; other types of values are stored in
148 * gmx::SelectionTreeElement::v.
149 * Constants that appear as parameters for selection methods are not present
150 * in the selection tree unless they have \ref GROUP_VALUE.
151 * \ref SEL_CONST elements have no children.
154 * \subsection selparser_tree_method Method evaluation elements
156 * \ref SEL_EXPRESSION and \ref SEL_MODIFIER elements are treated very
157 * similarly. The \c gmx_ana_selmethod_t structure corresponding to the
158 * evaluation method is in gmx::SelectionTreeElement::method, and the method
159 * data in gmx::SelectionTreeElement::mdata has been allocated using
160 * sel_datafunc().
161 * If a non-standard reference position type was set,
162 * gmx::SelectionTreeElement::pc has also been created, but only the type has
163 * been set.
164 * All children of these elements are of the type \ref SEL_SUBEXPRREF, and
165 * each describes a selection that needs to be evaluated to obtain a value
166 * for one parameter of the method.
167 * No children are present for parameters that were given a constant
168 * non-\ref GROUP_VALUE value.
169 * The children are sorted in the order in which the parameters appear in the
170 * \ref gmx_ana_selmethod_t structure.
172 * In addition to actual selection keywords, \ref SEL_EXPRESSION elements
173 * are used internally to implement numerical comparisons (e.g., "x < 5")
174 * and keyword matching (e.g., "resnr 1 to 3" or "name CA").
177 * \subsection selparser_tree_subexpr Subexpression elements
179 * \ref SEL_SUBEXPR elements only appear for variables, as described above.
180 * gmx::SelectionTreeElement::name points to the name of the variable (from the
181 * \ref SEL_ROOT element).
182 * The element always has exactly one child, which represents the value of
183 * the variable.
185 * \ref SEL_SUBEXPRREF elements are used for two purposes:
186 * - Variable references that need to be evaluated (i.e., there is a
187 * \ref SEL_SUBEXPR element for the variable) are represented using
188 * \ref SEL_SUBEXPRREF elements.
189 * In this case, gmx::SelectionTreeElement::param is NULL, and the first and
190 * only child of the element is the \ref SEL_SUBEXPR element of the
191 * variable.
192 * Such references can appear anywhere where the variable value
193 * (the child of the \ref SEL_SUBEXPR element) would be valid.
194 * - Children of \ref SEL_EXPRESSION and \ref SEL_MODIFIER elements are
195 * always of this type. For these elements, gmx::SelectionTreeElement::param
196 * is initialized to point to the parameter that receives the value from
197 * the expression.
198 * Each such element has exactly one child, which can be of any type;
199 * the \ref SEL_SUBEXPR element of a variable is used if the value comes
200 * from a variable, otherwise the child type is not \ref SEL_SUBEXPR.
203 * \subsection selparser_tree_bool Boolean elements
205 * One \ref SEL_BOOLEAN element is created for each boolean keyword in the
206 * input, and the tree structure represents the evaluation order.
207 * The gmx::SelectionTreeElement::boolt type gives the type of the operation.
208 * Each element has exactly two children (one for \ref BOOL_NOT elements),
209 * which are in the order given in the input.
210 * The children always have \ref GROUP_VALUE, but different element types
211 * are possible.
214 * \subsection selparser_tree_arith Arithmetic elements
216 * One \ref SEL_ARITHMETIC element is created for each arithmetic operation in
217 * the input, and the tree structure represents the evaluation order.
218 * The gmx::SelectionTreeElement::optype type gives the name of the operation.
219 * Each element has exactly two children (one for unary negation elements),
220 * which are in the order given in the input.
222 #include <stdio.h>
223 #include <stdarg.h>
225 #include <boost/exception_ptr.hpp>
226 #include <boost/shared_ptr.hpp>
228 #include "gromacs/fileio/futil.h"
229 #include "gromacs/selection/poscalc.h"
230 #include "gromacs/selection/selection.h"
231 #include "gromacs/selection/selmethod.h"
232 #include "gromacs/utility/exceptions.h"
233 #include "gromacs/utility/file.h"
234 #include "gromacs/utility/messagestringcollector.h"
235 #include "gromacs/utility/smalloc.h"
236 #include "gromacs/utility/stringutil.h"
238 #include "keywords.h"
239 #include "parsetree.h"
240 #include "selectioncollection-impl.h"
241 #include "selelem.h"
242 #include "symrec.h"
244 #include "scanner.h"
246 using gmx::SelectionParserValue;
247 using gmx::SelectionParserValueList;
248 using gmx::SelectionParserValueListPointer;
249 using gmx::SelectionParserParameter;
250 using gmx::SelectionParserParameterList;
251 using gmx::SelectionParserParameterListPointer;
252 using gmx::SelectionParserValue;
253 using gmx::SelectionTreeElement;
254 using gmx::SelectionTreeElementPointer;
256 void
257 _gmx_selparser_error(yyscan_t scanner, const char *fmt, ...)
259 gmx::MessageStringCollector *errors = _gmx_sel_lexer_error_reporter(scanner);
260 // FIXME: Use an arbitrary length buffer.
261 char buf[1024];
262 va_list ap;
263 va_start(ap, fmt);
264 vsnprintf(buf, 1024, fmt, ap);
265 va_end(ap);
266 errors->append(buf);
269 bool
270 _gmx_selparser_handle_exception(yyscan_t scanner, const std::exception &ex)
272 if (dynamic_cast<const gmx::UserInputError *>(&ex) != NULL)
274 // TODO: Consider whether also the non-interactive parser should
275 // postpone the exception such that the whole selection can be added as
276 // context.
277 if (_gmx_sel_is_lexer_interactive(scanner))
279 // TODO: Handle exceptions that printing the message may produce.
280 gmx::formatExceptionMessageToFile(stderr, ex);
281 return true;
284 _gmx_sel_lexer_set_exception(scanner, boost::current_exception());
285 return false;
288 namespace gmx
291 /********************************************************************
292 * SelectionParserValue
295 SelectionParserValue::SelectionParserValue(e_selvalue_t type)
296 : type(type)
298 memset(&u, 0, sizeof(u));
301 SelectionParserValue::SelectionParserValue(
302 const SelectionTreeElementPointer &expr)
303 : type(expr->v.type), expr(expr)
305 memset(&u, 0, sizeof(u));
308 /********************************************************************
309 * SelectionParserParameter
312 SelectionParserParameter::SelectionParserParameter(
313 const char *name,
314 SelectionParserValueListPointer values)
315 : name_(name != NULL ? name : ""),
316 values_(values ? move(values)
317 : SelectionParserValueListPointer(new SelectionParserValueList))
321 } // namespace gmx
324 * \param[in,out] sel Root of the selection element tree to initialize.
326 * Propagates the \ref SEL_DYNAMIC flag from the children of \p sel to \p sel
327 * (if any child of \p sel is dynamic, \p sel is also marked as such).
328 * The \ref SEL_DYNAMIC flag is also set for \ref SEL_EXPRESSION elements with
329 * a dynamic method.
330 * Also, sets one of the \ref SEL_SINGLEVAL, \ref SEL_ATOMVAL, or
331 * \ref SEL_VARNUMVAL flags, either based on the children or on the type of
332 * the selection method.
333 * If the types of the children conflict, an error is returned.
335 * The flags of the children of \p sel are also updated if not done earlier.
336 * The flags are initialized only once for any element; if \ref SEL_FLAGSSET
337 * is set for an element, the function returns immediately, and the recursive
338 * operation does not descend beyond such elements.
340 void
341 _gmx_selelem_update_flags(const gmx::SelectionTreeElementPointer &sel)
343 bool bUseChildType = false;
344 bool bOnlySingleChildren;
346 /* Return if the flags have already been set */
347 if (sel->flags & SEL_FLAGSSET)
349 return;
351 /* Set the flags based on the current element type */
352 switch (sel->type)
354 case SEL_CONST:
355 case SEL_GROUPREF:
356 sel->flags |= SEL_SINGLEVAL;
357 bUseChildType = false;
358 break;
360 case SEL_EXPRESSION:
361 if (sel->u.expr.method->flags & SMETH_DYNAMIC)
363 sel->flags |= SEL_DYNAMIC;
365 if (sel->u.expr.method->flags & SMETH_SINGLEVAL)
367 sel->flags |= SEL_SINGLEVAL;
369 else if (sel->u.expr.method->flags & SMETH_VARNUMVAL)
371 sel->flags |= SEL_VARNUMVAL;
373 else
375 sel->flags |= SEL_ATOMVAL;
377 bUseChildType = false;
378 break;
380 case SEL_ARITHMETIC:
381 sel->flags |= SEL_ATOMVAL;
382 bUseChildType = false;
383 break;
385 case SEL_MODIFIER:
386 if (sel->v.type != NO_VALUE)
388 sel->flags |= SEL_VARNUMVAL;
390 bUseChildType = false;
391 break;
393 case SEL_ROOT:
394 bUseChildType = false;
395 break;
397 case SEL_BOOLEAN:
398 case SEL_SUBEXPR:
399 case SEL_SUBEXPRREF:
400 bUseChildType = true;
401 break;
403 /* Loop through children to propagate their flags upwards */
404 bOnlySingleChildren = true;
405 SelectionTreeElementPointer child = sel->child;
406 while (child)
408 /* Update the child */
409 _gmx_selelem_update_flags(child);
410 /* Propagate the dynamic and unsorted flags */
411 sel->flags |= (child->flags & (SEL_DYNAMIC | SEL_UNSORTED));
412 /* Propagate the type flag if necessary and check for problems */
413 if (bUseChildType)
415 if ((sel->flags & SEL_VALTYPEMASK)
416 && !(sel->flags & child->flags & SEL_VALTYPEMASK))
418 // TODO: Recollect when this is triggered, and whether the type
419 // is appropriate.
420 GMX_THROW(gmx::InvalidInputError("Invalid combination of selection expressions"));
422 sel->flags |= (child->flags & SEL_VALTYPEMASK);
424 if (!(child->flags & SEL_SINGLEVAL))
426 bOnlySingleChildren = false;
429 child = child->next;
431 /* For arithmetic expressions consisting only of single values,
432 * the result is also a single value. */
433 if (sel->type == SEL_ARITHMETIC && bOnlySingleChildren)
435 sel->flags = (sel->flags & ~SEL_VALTYPEMASK) | SEL_SINGLEVAL;
437 /* For root elements, the type should be propagated here, after the
438 * children have been updated. */
439 if (sel->type == SEL_ROOT)
441 GMX_ASSERT(sel->child, "Root elements should always have a child");
442 sel->flags |= (sel->child->flags & SEL_VALTYPEMASK);
444 /* Mark that the flags are set */
445 sel->flags |= SEL_FLAGSSET;
449 * \param[in,out] sel Selection element to initialize.
450 * \param[in] scanner Scanner data structure.
452 * A deep copy of the parameters is made to allow several
453 * expressions with the same method to coexist peacefully.
454 * Calls sel_datafunc() if one is specified for the method.
456 void
457 _gmx_selelem_init_method_params(const gmx::SelectionTreeElementPointer &sel,
458 yyscan_t scanner)
460 int nparams;
461 gmx_ana_selparam_t *orgparam;
462 gmx_ana_selparam_t *param;
463 int i;
464 void *mdata;
466 nparams = sel->u.expr.method->nparams;
467 orgparam = sel->u.expr.method->param;
468 snew(param, nparams);
469 memcpy(param, orgparam, nparams*sizeof(gmx_ana_selparam_t));
470 for (i = 0; i < nparams; ++i)
472 param[i].flags &= ~SPAR_SET;
473 _gmx_selvalue_setstore(&param[i].val, NULL);
474 if (param[i].flags & SPAR_VARNUM)
476 param[i].val.nr = -1;
478 /* Duplicate the enum value array if it is given statically */
479 if ((param[i].flags & SPAR_ENUMVAL) && orgparam[i].val.u.ptr != NULL)
481 int n;
483 /* Count the values */
484 n = 1;
485 while (orgparam[i].val.u.s[n] != NULL)
487 ++n;
489 _gmx_selvalue_reserve(&param[i].val, n+1);
490 memcpy(param[i].val.u.s, orgparam[i].val.u.s,
491 (n+1)*sizeof(param[i].val.u.s[0]));
494 mdata = NULL;
495 if (sel->u.expr.method->init_data)
497 mdata = sel->u.expr.method->init_data(nparams, param);
499 if (sel->u.expr.method->set_poscoll)
501 gmx_ana_selcollection_t *sc = _gmx_sel_lexer_selcollection(scanner);
503 sel->u.expr.method->set_poscoll(&sc->pcc, mdata);
505 /* Store the values */
506 sel->u.expr.method->param = param;
507 sel->u.expr.mdata = mdata;
511 * \param[in,out] sel Selection element to initialize.
512 * \param[in] method Selection method to set.
513 * \param[in] scanner Scanner data structure.
515 * Makes a copy of \p method and stores it in \p sel->u.expr.method,
516 * and calls _gmx_selelem_init_method_params();
518 void
519 _gmx_selelem_set_method(const gmx::SelectionTreeElementPointer &sel,
520 gmx_ana_selmethod_t *method,
521 yyscan_t scanner)
523 _gmx_selelem_set_vtype(sel, method->type);
524 sel->setName(method->name);
525 snew(sel->u.expr.method, 1);
526 memcpy(sel->u.expr.method, method, sizeof(gmx_ana_selmethod_t));
527 _gmx_selelem_init_method_params(sel, scanner);
530 /*! \brief
531 * Initializes the reference position calculation for a \ref SEL_EXPRESSION
532 * element.
534 * \param[in,out] pcc Position calculation collection to use.
535 * \param[in,out] sel Selection element to initialize.
536 * \param[in] rpost Reference position type to use (NULL = default).
537 * \param[in] scanner Scanner data structure.
538 * \returns 0 on success, a non-zero error code on error.
540 static void
541 set_refpos_type(gmx::PositionCalculationCollection *pcc,
542 const SelectionTreeElementPointer &sel,
543 const char *rpost, yyscan_t scanner)
545 if (!rpost)
547 return;
550 if (sel->u.expr.method->pupdate)
552 /* By default, use whole residues/molecules. */
553 sel->u.expr.pc
554 = pcc->createCalculationFromEnum(rpost, POS_COMPLWHOLE);
556 else
558 // TODO: Should this be treated as a real error?
559 _gmx_selparser_error(scanner, "modifier '%s' is not applicable for '%s'",
560 rpost, sel->u.expr.method->name);
564 gmx::SelectionTreeElementPointer
565 _gmx_sel_init_arithmetic(const gmx::SelectionTreeElementPointer &left,
566 const gmx::SelectionTreeElementPointer &right,
567 char op, yyscan_t /* scanner */)
569 SelectionTreeElementPointer sel(new SelectionTreeElement(SEL_ARITHMETIC));
570 sel->v.type = REAL_VALUE;
571 switch (op)
573 case '+': sel->u.arith.type = ARITH_PLUS; break;
574 case '-': sel->u.arith.type = (right ? ARITH_MINUS : ARITH_NEG); break;
575 case '*': sel->u.arith.type = ARITH_MULT; break;
576 case '/': sel->u.arith.type = ARITH_DIV; break;
577 case '^': sel->u.arith.type = ARITH_EXP; break;
579 char buf[2];
580 buf[0] = op;
581 buf[1] = 0;
582 sel->setName(buf);
583 sel->u.arith.opstr = strdup(buf);
584 sel->child = left;
585 sel->child->next = right;
586 return sel;
590 * \param[in] left Selection element for the left hand side.
591 * \param[in] right Selection element for the right hand side.
592 * \param[in] cmpop String representation of the comparison operator.
593 * \param[in] scanner Scanner data structure.
594 * \returns The created selection element.
596 * This function handles the creation of a gmx::SelectionTreeElement object for
597 * comparison expressions.
599 SelectionTreeElementPointer
600 _gmx_sel_init_comparison(const gmx::SelectionTreeElementPointer &left,
601 const gmx::SelectionTreeElementPointer &right,
602 const char *cmpop, yyscan_t scanner)
604 gmx::MessageStringCollector *errors = _gmx_sel_lexer_error_reporter(scanner);
605 gmx::MessageStringContext context(errors, "In comparison initialization");
607 SelectionTreeElementPointer sel(new SelectionTreeElement(SEL_EXPRESSION));
608 _gmx_selelem_set_method(sel, &sm_compare, scanner);
610 SelectionParserParameterList params;
611 const char *name;
612 // Create the parameter for the left expression.
613 name = left->v.type == INT_VALUE ? "int1" : "real1";
614 params.push_back(SelectionParserParameter::createFromExpression(name, left));
615 // Create the parameter for the right expression.
616 name = right->v.type == INT_VALUE ? "int2" : "real2";
617 params.push_back(SelectionParserParameter::createFromExpression(name, right));
618 // Create the parameter for the operator.
619 params.push_back(
620 SelectionParserParameter::create(
621 "op", SelectionParserValue::createString(cmpop)));
622 if (!_gmx_sel_parse_params(params, sel->u.expr.method->nparams,
623 sel->u.expr.method->param, sel, scanner))
625 return SelectionTreeElementPointer();
628 return sel;
631 /*! \brief
632 * Implementation method for keyword expression creation.
634 * \param[in] method Method to use.
635 * \param[in] matchType String matching type (only used if \p method is
636 * a string keyword and \p args is not empty.
637 * \param[in] args Pointer to the first argument.
638 * \param[in] rpost Reference position type to use (NULL = default).
639 * \param[in] scanner Scanner data structure.
640 * \returns The created selection element.
642 * This function handles the creation of a gmx::SelectionTreeElement object for
643 * selection methods that do not take parameters.
645 static SelectionTreeElementPointer
646 init_keyword_internal(gmx_ana_selmethod_t *method,
647 gmx::SelectionStringMatchType matchType,
648 SelectionParserValueListPointer args,
649 const char *rpost, yyscan_t scanner)
651 gmx_ana_selcollection_t *sc = _gmx_sel_lexer_selcollection(scanner);
653 gmx::MessageStringCollector *errors = _gmx_sel_lexer_error_reporter(scanner);
654 char buf[128];
655 sprintf(buf, "In keyword '%s'", method->name);
656 gmx::MessageStringContext context(errors, buf);
658 if (method->nparams > 0)
660 // TODO: Would assert be better?
661 GMX_THROW(gmx::InternalError(
662 "Keyword initialization called with non-keyword method"));
665 SelectionTreeElementPointer root(new SelectionTreeElement(SEL_EXPRESSION));
666 SelectionTreeElementPointer child = root;
667 _gmx_selelem_set_method(child, method, scanner);
669 /* Initialize the evaluation of keyword matching if values are provided */
670 if (args)
672 gmx_ana_selmethod_t *kwmethod;
673 switch (method->type)
675 case INT_VALUE: kwmethod = &sm_keyword_int; break;
676 case REAL_VALUE: kwmethod = &sm_keyword_real; break;
677 case STR_VALUE: kwmethod = &sm_keyword_str; break;
678 default:
679 GMX_THROW(gmx::InternalError(
680 "Unknown type for keyword selection"));
682 /* Initialize the selection element */
683 root.reset(new SelectionTreeElement(SEL_EXPRESSION));
684 _gmx_selelem_set_method(root, kwmethod, scanner);
685 if (method->type == STR_VALUE)
687 _gmx_selelem_set_kwstr_match_type(root, matchType);
689 SelectionParserParameterList params;
690 params.push_back(
691 SelectionParserParameter::createFromExpression(NULL, child));
692 params.push_back(SelectionParserParameter::create(NULL, move(args)));
693 if (!_gmx_sel_parse_params(params, root->u.expr.method->nparams,
694 root->u.expr.method->param, root, scanner))
696 return SelectionTreeElementPointer();
699 set_refpos_type(&sc->pcc, child, rpost, scanner);
701 return root;
705 * \param[in] method Method to use.
706 * \param[in] args Pointer to the first argument.
707 * \param[in] rpost Reference position type to use (NULL = default).
708 * \param[in] scanner Scanner data structure.
709 * \returns The created selection element.
711 * This function handles the creation of a gmx::SelectionTreeElement object for
712 * selection methods that do not take parameters.
714 SelectionTreeElementPointer
715 _gmx_sel_init_keyword(gmx_ana_selmethod_t *method,
716 gmx::SelectionParserValueListPointer args,
717 const char *rpost, yyscan_t scanner)
719 return init_keyword_internal(method, gmx::eStringMatchType_Auto, move(args),
720 rpost, scanner);
724 * \param[in] method Method to use.
725 * \param[in] matchType String matching type.
726 * \param[in] args Pointer to the first argument.
727 * \param[in] rpost Reference position type to use (NULL = default).
728 * \param[in] scanner Scanner data structure.
729 * \returns The created selection element.
731 * This function handles the creation of a gmx::SelectionTreeElement object for
732 * keyword string matching.
734 SelectionTreeElementPointer
735 _gmx_sel_init_keyword_strmatch(gmx_ana_selmethod_t *method,
736 gmx::SelectionStringMatchType matchType,
737 gmx::SelectionParserValueListPointer args,
738 const char *rpost, yyscan_t scanner)
740 GMX_RELEASE_ASSERT(method->type == STR_VALUE,
741 "String keyword method called for a non-string-valued method");
742 GMX_RELEASE_ASSERT(args && !args->empty(),
743 "String keyword matching method called without any values");
744 return init_keyword_internal(method, matchType, move(args), rpost, scanner);
748 * \param[in] method Method to use for initialization.
749 * \param[in] params Pointer to the first parameter.
750 * \param[in] rpost Reference position type to use (NULL = default).
751 * \param[in] scanner Scanner data structure.
752 * \returns The created selection element.
754 * This function handles the creation of a gmx::SelectionTreeElement object for
755 * selection methods that take parameters.
757 * Part of the behavior of the \c same selection keyword is hardcoded into
758 * this function (or rather, into _gmx_selelem_custom_init_same()) to allow the
759 * use of any keyword in \c "same KEYWORD as" without requiring special
760 * handling somewhere else (or sacrificing the simple syntax).
762 SelectionTreeElementPointer
763 _gmx_sel_init_method(gmx_ana_selmethod_t *method,
764 gmx::SelectionParserParameterListPointer params,
765 const char *rpost, yyscan_t scanner)
767 gmx_ana_selcollection_t *sc = _gmx_sel_lexer_selcollection(scanner);
768 int rc;
770 gmx::MessageStringCollector *errors = _gmx_sel_lexer_error_reporter(scanner);
771 char buf[128];
772 sprintf(buf, "In keyword '%s'", method->name);
773 gmx::MessageStringContext context(errors, buf);
775 _gmx_sel_finish_method(scanner);
776 /* The "same" keyword needs some custom massaging of the parameters. */
777 rc = _gmx_selelem_custom_init_same(&method, params, scanner);
778 if (rc != 0)
780 return SelectionTreeElementPointer();
782 SelectionTreeElementPointer root(new SelectionTreeElement(SEL_EXPRESSION));
783 _gmx_selelem_set_method(root, method, scanner);
784 /* Process the parameters */
785 if (!_gmx_sel_parse_params(*params, root->u.expr.method->nparams,
786 root->u.expr.method->param, root, scanner))
788 return SelectionTreeElementPointer();
790 set_refpos_type(&sc->pcc, root, rpost, scanner);
792 return root;
796 * \param[in] method Modifier to use for initialization.
797 * \param[in] params Pointer to the first parameter.
798 * \param[in] sel Selection element that the modifier should act on.
799 * \param[in] scanner Scanner data structure.
800 * \returns The created selection element.
802 * This function handles the creation of a gmx::SelectionTreeElement object for
803 * selection modifiers.
805 SelectionTreeElementPointer
806 _gmx_sel_init_modifier(gmx_ana_selmethod_t *method,
807 gmx::SelectionParserParameterListPointer params,
808 const gmx::SelectionTreeElementPointer &sel,
809 yyscan_t scanner)
811 gmx::MessageStringCollector *errors = _gmx_sel_lexer_error_reporter(scanner);
812 char buf[128];
813 sprintf(buf, "In keyword '%s'", method->name);
814 gmx::MessageStringContext context(errors, buf);
816 _gmx_sel_finish_method(scanner);
817 SelectionTreeElementPointer modifier(new SelectionTreeElement(SEL_MODIFIER));
818 _gmx_selelem_set_method(modifier, method, scanner);
819 SelectionTreeElementPointer root;
820 if (method->type == NO_VALUE)
822 SelectionTreeElementPointer child = sel;
823 while (child->next)
825 child = child->next;
827 child->next = modifier;
828 root = sel;
830 else
832 params->push_front(
833 SelectionParserParameter::createFromExpression(NULL, sel));
834 root = modifier;
836 /* Process the parameters */
837 if (!_gmx_sel_parse_params(*params, modifier->u.expr.method->nparams,
838 modifier->u.expr.method->param, modifier, scanner))
840 return SelectionTreeElementPointer();
843 return root;
847 * \param[in] expr Input selection element for the position calculation.
848 * \param[in] type Reference position type or NULL for default.
849 * \param[in] scanner Scanner data structure.
850 * \returns The created selection element.
852 * This function handles the creation of a gmx::SelectionTreeElement object for
853 * evaluation of reference positions.
855 SelectionTreeElementPointer
856 _gmx_sel_init_position(const gmx::SelectionTreeElementPointer &expr,
857 const char *type, yyscan_t scanner)
859 gmx::MessageStringCollector *errors = _gmx_sel_lexer_error_reporter(scanner);
860 char buf[128];
861 sprintf(buf, "In position evaluation");
862 gmx::MessageStringContext context(errors, buf);
864 SelectionTreeElementPointer root(new SelectionTreeElement(SEL_EXPRESSION));
865 _gmx_selelem_set_method(root, &sm_keyword_pos, scanner);
866 _gmx_selelem_set_kwpos_type(root.get(), type);
867 /* Create the parameters for the parameter parser. */
868 SelectionParserParameterList params;
869 params.push_back(SelectionParserParameter::createFromExpression(NULL, expr));
870 /* Parse the parameters. */
871 if (!_gmx_sel_parse_params(params, root->u.expr.method->nparams,
872 root->u.expr.method->param, root, scanner))
874 return SelectionTreeElementPointer();
877 return root;
881 * \param[in] x,y,z Coordinates for the position.
882 * \returns The creates selection element.
884 SelectionTreeElementPointer
885 _gmx_sel_init_const_position(real x, real y, real z)
887 rvec pos;
889 SelectionTreeElementPointer sel(new SelectionTreeElement(SEL_CONST));
890 _gmx_selelem_set_vtype(sel, POS_VALUE);
891 _gmx_selvalue_reserve(&sel->v, 1);
892 pos[XX] = x;
893 pos[YY] = y;
894 pos[ZZ] = z;
895 gmx_ana_pos_init_const(sel->v.u.p, pos);
896 return sel;
900 * \param[in] name Name of an index group to search for.
901 * \param[in] scanner Scanner data structure.
902 * \returns The created selection element.
904 * See gmx_ana_indexgrps_find() for information on how \p name is matched
905 * against the index groups.
907 SelectionTreeElementPointer
908 _gmx_sel_init_group_by_name(const char *name, yyscan_t scanner)
911 SelectionTreeElementPointer sel(new SelectionTreeElement(SEL_GROUPREF));
912 _gmx_selelem_set_vtype(sel, GROUP_VALUE);
913 sel->setName(gmx::formatString("group \"%s\"", name));
914 sel->u.gref.name = strdup(name);
915 sel->u.gref.id = -1;
917 if (_gmx_sel_lexer_has_groups_set(scanner))
919 gmx_ana_indexgrps_t *grps = _gmx_sel_lexer_indexgrps(scanner);
920 sel->resolveIndexGroupReference(grps);
923 return sel;
927 * \param[in] id Zero-based index number of the group to extract.
928 * \param[in] scanner Scanner data structure.
929 * \returns The created selection element.
931 SelectionTreeElementPointer
932 _gmx_sel_init_group_by_id(int id, yyscan_t scanner)
934 SelectionTreeElementPointer sel(new SelectionTreeElement(SEL_GROUPREF));
935 _gmx_selelem_set_vtype(sel, GROUP_VALUE);
936 sel->setName(gmx::formatString("group %d", id));
937 sel->u.gref.name = NULL;
938 sel->u.gref.id = id;
940 if (_gmx_sel_lexer_has_groups_set(scanner))
942 gmx_ana_indexgrps_t *grps = _gmx_sel_lexer_indexgrps(scanner);
943 sel->resolveIndexGroupReference(grps);
946 return sel;
950 * \param[in,out] sel Value of the variable.
951 * \returns The created selection element that references \p sel.
953 * The reference count of \p sel is updated, but no other modifications are
954 * made.
956 SelectionTreeElementPointer
957 _gmx_sel_init_variable_ref(const gmx::SelectionTreeElementPointer &sel)
959 SelectionTreeElementPointer ref;
961 if (sel->v.type == POS_VALUE && sel->type == SEL_CONST)
963 ref = sel;
965 else
967 ref.reset(new SelectionTreeElement(SEL_SUBEXPRREF));
968 _gmx_selelem_set_vtype(ref, sel->v.type);
969 ref->setName(sel->name());
970 ref->child = sel;
972 return ref;
976 * \param[in] name Name for the selection
977 * (if NULL, a default name is constructed).
978 * \param[in] sel The selection element that evaluates the selection.
979 * \param scanner Scanner data structure.
980 * \returns The created root selection element.
982 * This function handles the creation of root (\ref SEL_ROOT)
983 * gmx::SelectionTreeElement objects for selections.
985 SelectionTreeElementPointer
986 _gmx_sel_init_selection(const char *name,
987 const gmx::SelectionTreeElementPointer &sel,
988 yyscan_t scanner)
990 if (sel->v.type != POS_VALUE)
992 /* FIXME: Better handling of this error */
993 GMX_THROW(gmx::InternalError(
994 "Each selection must evaluate to a position"));
997 SelectionTreeElementPointer root(new SelectionTreeElement(SEL_ROOT));
998 root->child = sel;
999 if (name)
1001 root->setName(name);
1003 /* Update the flags */
1004 _gmx_selelem_update_flags(root);
1006 root->fillNameIfMissing(_gmx_sel_lexer_pselstr(scanner));
1008 /* Print out some information if the parser is interactive */
1009 if (_gmx_sel_is_lexer_interactive(scanner))
1011 fprintf(stderr, "Selection '%s' parsed\n",
1012 _gmx_sel_lexer_pselstr(scanner));
1015 return root;
1020 * \param[in] name Name of the variable.
1021 * \param[in] expr The selection element that evaluates the variable.
1022 * \param scanner Scanner data structure.
1023 * \returns The created root selection element.
1025 * This function handles the creation of root gmx::SelectionTreeElement objects
1026 * for variable assignments. A \ref SEL_ROOT element and a \ref SEL_SUBEXPR
1027 * element are both created.
1029 SelectionTreeElementPointer
1030 _gmx_sel_assign_variable(const char *name,
1031 const gmx::SelectionTreeElementPointer &expr,
1032 yyscan_t scanner)
1034 gmx_ana_selcollection_t *sc = _gmx_sel_lexer_selcollection(scanner);
1035 const char *pselstr = _gmx_sel_lexer_pselstr(scanner);
1036 SelectionTreeElementPointer root;
1038 _gmx_selelem_update_flags(expr);
1039 /* Check if this is a constant non-group value */
1040 if (expr->type == SEL_CONST && expr->v.type != GROUP_VALUE)
1042 /* If so, just assign the constant value to the variable */
1043 sc->symtab->addVariable(name, expr);
1044 goto finish;
1046 /* Check if we are assigning a variable to another variable */
1047 if (expr->type == SEL_SUBEXPRREF)
1049 /* If so, make a simple alias */
1050 sc->symtab->addVariable(name, expr->child);
1051 goto finish;
1053 /* Create the root element */
1054 root.reset(new SelectionTreeElement(SEL_ROOT));
1055 root->setName(name);
1056 /* Create the subexpression element */
1057 root->child.reset(new SelectionTreeElement(SEL_SUBEXPR));
1058 root->child->setName(name);
1059 _gmx_selelem_set_vtype(root->child, expr->v.type);
1060 root->child->child = expr;
1061 /* Update flags */
1062 _gmx_selelem_update_flags(root);
1063 /* Add the variable to the symbol table */
1064 sc->symtab->addVariable(name, root->child);
1065 finish:
1066 srenew(sc->varstrs, sc->nvars + 1);
1067 sc->varstrs[sc->nvars] = strdup(pselstr);
1068 ++sc->nvars;
1069 if (_gmx_sel_is_lexer_interactive(scanner))
1071 fprintf(stderr, "Variable '%s' parsed\n", pselstr);
1073 return root;
1077 * \param sel Selection to append (can be NULL, in which
1078 * case nothing is done).
1079 * \param last Last selection, or NULL if not present or not known.
1080 * \param scanner Scanner data structure.
1081 * \returns The last selection after the append.
1083 * Appends \p sel after the last root element, and returns either \p sel
1084 * (if it was non-NULL) or the last element (if \p sel was NULL).
1086 SelectionTreeElementPointer
1087 _gmx_sel_append_selection(const gmx::SelectionTreeElementPointer &sel,
1088 gmx::SelectionTreeElementPointer last,
1089 yyscan_t scanner)
1091 gmx_ana_selcollection_t *sc = _gmx_sel_lexer_selcollection(scanner);
1093 /* Append sel after last, or the last element of sc if last is NULL */
1094 if (last)
1096 last->next = sel;
1098 else
1100 if (sc->root)
1102 last = sc->root;
1103 while (last->next)
1105 last = last->next;
1107 last->next = sel;
1109 else
1111 sc->root = sel;
1114 /* Initialize a selection object if necessary */
1115 if (sel)
1117 last = sel;
1118 /* Add the new selection to the collection if it is not a variable. */
1119 if (sel->child->type != SEL_SUBEXPR)
1121 gmx::SelectionDataPointer selPtr(
1122 new gmx::internal::SelectionData(
1123 sel.get(), _gmx_sel_lexer_pselstr(scanner)));
1124 sc->sel.push_back(gmx::move(selPtr));
1127 /* Clear the selection string now that we've saved it */
1128 _gmx_sel_lexer_clear_pselstr(scanner);
1129 return last;
1133 * \param[in] scanner Scanner data structure.
1134 * \returns true if the parser should finish, false if parsing should
1135 * continue.
1137 * This function is called always after _gmx_sel_append_selection() to
1138 * check whether a sufficient number of selections has already been provided.
1139 * This is used to terminate interactive parsers when the correct number of
1140 * selections has been provided.
1142 bool
1143 _gmx_sel_parser_should_finish(yyscan_t scanner)
1145 gmx_ana_selcollection_t *sc = _gmx_sel_lexer_selcollection(scanner);
1146 return (int)sc->sel.size() == _gmx_sel_lexer_exp_selcount(scanner);