Merge release-4-6 into release-5-0
[alexxy/gromacs.git] / src / gromacs / selection / parsetree.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 /*! \internal \file
36  * \brief
37  * Implements functions in parsetree.h.
38  *
39  * \author Teemu Murtola <teemu.murtola@gmail.com>
40  * \ingroup module_selection
41  */
42 /*! \internal
43  * \page page_module_selection_parser Selection parsing
44  *
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().
72  *
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.
87  *
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
93  *
94  *
95  * \section selparser_tree Element tree constructed by the parser
96  *
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.
106  *
107  *
108  * \subsection selparser_tree_root Root elements
109  *
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.
116  *  .
117  * The \ref SEL_ROOT elements are linked together in a chain in the same order
118  * as in the input.
119  *
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.
129  *  .
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.
135  *
136  *
137  * \subsection selparser_tree_const Constant elements
138  *
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,
145  *  .
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.
152  *
153  *
154  * \subsection selparser_tree_method Method evaluation elements
155  *
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.
171  *
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").
175  *
176  *
177  * \subsection selparser_tree_subexpr Subexpression elements
178  *
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.
184  *
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.
201  *
202  *
203  * \subsection selparser_tree_bool Boolean elements
204  *
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.
212  *
213  *
214  * \subsection selparser_tree_arith Arithmetic elements
215  *
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.
221  */
222 #include <stdio.h>
223 #include <stdarg.h>
224
225 #include <boost/exception_ptr.hpp>
226 #include <boost/shared_ptr.hpp>
227
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"
237
238 #include "keywords.h"
239 #include "parsetree.h"
240 #include "selectioncollection-impl.h"
241 #include "selelem.h"
242 #include "symrec.h"
243
244 #include "scanner.h"
245
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;
255
256 void
257 _gmx_selparser_error(yyscan_t scanner, const char *fmt, ...)
258 {
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     vsprintf(buf, fmt, ap);
265     va_end(ap);
266     errors->append(buf);
267 }
268
269 bool
270 _gmx_selparser_handle_exception(yyscan_t scanner, const std::exception &ex)
271 {
272     if (dynamic_cast<const gmx::UserInputError *>(&ex) != NULL)
273     {
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))
278         {
279             // TODO: Handle exceptions that printing the message may produce.
280             gmx::formatExceptionMessageToFile(stderr, ex);
281             return true;
282         }
283     }
284     _gmx_sel_lexer_set_exception(scanner, boost::current_exception());
285     return false;
286 }
287
288 namespace gmx
289 {
290
291 /********************************************************************
292  * SelectionParserValue
293  */
294
295 SelectionParserValue::SelectionParserValue(e_selvalue_t type)
296     : type(type)
297 {
298     memset(&u, 0, sizeof(u));
299 }
300
301 SelectionParserValue::SelectionParserValue(
302         const SelectionTreeElementPointer &expr)
303     : type(expr->v.type), expr(expr)
304 {
305     memset(&u, 0, sizeof(u));
306 }
307
308 /********************************************************************
309  * SelectionParserParameter
310  */
311
312 SelectionParserParameter::SelectionParserParameter(
313         const char                     *name,
314         SelectionParserValueListPointer values)
315     : name_(name != NULL ? name : ""),
316       values_(values ? move(values)
317               : SelectionParserValueListPointer(new SelectionParserValueList))
318 {
319 }
320
321 } // namespace gmx
322
323 /*!
324  * \param[in,out] sel  Root of the selection element tree to initialize.
325  * \param[in]     scanner Scanner data structure.
326  * \returns       0 on success, an error code on error.
327  *
328  * Propagates the \ref SEL_DYNAMIC flag from the children of \p sel to \p sel
329  * (if any child of \p sel is dynamic, \p sel is also marked as such).
330  * The \ref SEL_DYNAMIC flag is also set for \ref SEL_EXPRESSION elements with
331  * a dynamic method.
332  * Also, sets one of the \ref SEL_SINGLEVAL, \ref SEL_ATOMVAL, or
333  * \ref SEL_VARNUMVAL flags, either based on the children or on the type of
334  * the selection method.
335  * If the types of the children conflict, an error is returned.
336  *
337  * The flags of the children of \p sel are also updated if not done earlier.
338  * The flags are initialized only once for any element; if \ref SEL_FLAGSSET
339  * is set for an element, the function returns immediately, and the recursive
340  * operation does not descend beyond such elements.
341  */
342 void
343 _gmx_selelem_update_flags(const gmx::SelectionTreeElementPointer &sel,
344                           yyscan_t                                scanner)
345 {
346     bool                bUseChildType = false;
347     bool                bOnlySingleChildren;
348
349     /* Return if the flags have already been set */
350     if (sel->flags & SEL_FLAGSSET)
351     {
352         return;
353     }
354     /* Set the flags based on the current element type */
355     switch (sel->type)
356     {
357         case SEL_CONST:
358         case SEL_GROUPREF:
359             sel->flags   |= SEL_SINGLEVAL;
360             bUseChildType = false;
361             break;
362
363         case SEL_EXPRESSION:
364             if (sel->u.expr.method->flags & SMETH_DYNAMIC)
365             {
366                 sel->flags |= SEL_DYNAMIC;
367             }
368             if (sel->u.expr.method->flags & SMETH_SINGLEVAL)
369             {
370                 sel->flags |= SEL_SINGLEVAL;
371             }
372             else if (sel->u.expr.method->flags & SMETH_VARNUMVAL)
373             {
374                 sel->flags |= SEL_VARNUMVAL;
375             }
376             else
377             {
378                 sel->flags |= SEL_ATOMVAL;
379             }
380             bUseChildType = false;
381             break;
382
383         case SEL_ARITHMETIC:
384             sel->flags   |= SEL_ATOMVAL;
385             bUseChildType = false;
386             break;
387
388         case SEL_MODIFIER:
389             if (sel->v.type != NO_VALUE)
390             {
391                 sel->flags |= SEL_VARNUMVAL;
392             }
393             bUseChildType = false;
394             break;
395
396         case SEL_ROOT:
397             bUseChildType = false;
398             break;
399
400         case SEL_BOOLEAN:
401         case SEL_SUBEXPR:
402         case SEL_SUBEXPRREF:
403             bUseChildType = true;
404             break;
405     }
406     /* Loop through children to propagate their flags upwards */
407     bOnlySingleChildren = true;
408     SelectionTreeElementPointer child = sel->child;
409     while (child)
410     {
411         /* Update the child */
412         _gmx_selelem_update_flags(child, scanner);
413         /* Propagate the dynamic flag */
414         sel->flags |= (child->flags & SEL_DYNAMIC);
415         /* Propagate the type flag if necessary and check for problems */
416         if (bUseChildType)
417         {
418             if ((sel->flags & SEL_VALTYPEMASK)
419                 && !(sel->flags & child->flags & SEL_VALTYPEMASK))
420             {
421                 _gmx_selparser_error(scanner, "invalid combination of selection expressions");
422                 // FIXME: Use an exception.
423                 return;
424             }
425             sel->flags |= (child->flags & SEL_VALTYPEMASK);
426         }
427         if (!(child->flags & SEL_SINGLEVAL))
428         {
429             bOnlySingleChildren = false;
430         }
431
432         child = child->next;
433     }
434     /* For arithmetic expressions consisting only of single values,
435      * the result is also a single value. */
436     if (sel->type == SEL_ARITHMETIC && bOnlySingleChildren)
437     {
438         sel->flags = (sel->flags & ~SEL_VALTYPEMASK) | SEL_SINGLEVAL;
439     }
440     /* For root elements, the type should be propagated here, after the
441      * children have been updated. */
442     if (sel->type == SEL_ROOT)
443     {
444         GMX_ASSERT(sel->child, "Root elements should always have a child");
445         sel->flags |= (sel->child->flags & SEL_VALTYPEMASK);
446     }
447     /* Mark that the flags are set */
448     sel->flags |= SEL_FLAGSSET;
449 }
450
451 /*!
452  * \param[in,out] sel    Selection element to initialize.
453  * \param[in]     scanner Scanner data structure.
454  *
455  * A deep copy of the parameters is made to allow several
456  * expressions with the same method to coexist peacefully.
457  * Calls sel_datafunc() if one is specified for the method.
458  */
459 void
460 _gmx_selelem_init_method_params(const gmx::SelectionTreeElementPointer &sel,
461                                 yyscan_t                                scanner)
462 {
463     int                 nparams;
464     gmx_ana_selparam_t *orgparam;
465     gmx_ana_selparam_t *param;
466     int                 i;
467     void               *mdata;
468
469     nparams   = sel->u.expr.method->nparams;
470     orgparam  = sel->u.expr.method->param;
471     snew(param, nparams);
472     memcpy(param, orgparam, nparams*sizeof(gmx_ana_selparam_t));
473     for (i = 0; i < nparams; ++i)
474     {
475         param[i].flags &= ~SPAR_SET;
476         _gmx_selvalue_setstore(&param[i].val, NULL);
477         if (param[i].flags & SPAR_VARNUM)
478         {
479             param[i].val.nr = -1;
480         }
481         /* Duplicate the enum value array if it is given statically */
482         if ((param[i].flags & SPAR_ENUMVAL) && orgparam[i].val.u.ptr != NULL)
483         {
484             int n;
485
486             /* Count the values */
487             n = 1;
488             while (orgparam[i].val.u.s[n] != NULL)
489             {
490                 ++n;
491             }
492             _gmx_selvalue_reserve(&param[i].val, n+1);
493             memcpy(param[i].val.u.s, orgparam[i].val.u.s,
494                    (n+1)*sizeof(param[i].val.u.s[0]));
495         }
496     }
497     mdata = NULL;
498     if (sel->u.expr.method->init_data)
499     {
500         mdata = sel->u.expr.method->init_data(nparams, param);
501     }
502     if (sel->u.expr.method->set_poscoll)
503     {
504         gmx_ana_selcollection_t *sc = _gmx_sel_lexer_selcollection(scanner);
505
506         sel->u.expr.method->set_poscoll(&sc->pcc, mdata);
507     }
508     /* Store the values */
509     sel->u.expr.method->param = param;
510     sel->u.expr.mdata         = mdata;
511 }
512
513 /*!
514  * \param[in,out] sel    Selection element to initialize.
515  * \param[in]     method Selection method to set.
516  * \param[in]     scanner Scanner data structure.
517  *
518  * Makes a copy of \p method and stores it in \p sel->u.expr.method,
519  * and calls _gmx_selelem_init_method_params();
520  */
521 void
522 _gmx_selelem_set_method(const gmx::SelectionTreeElementPointer &sel,
523                         gmx_ana_selmethod_t                    *method,
524                         yyscan_t                                scanner)
525 {
526     _gmx_selelem_set_vtype(sel, method->type);
527     sel->setName(method->name);
528     snew(sel->u.expr.method, 1);
529     memcpy(sel->u.expr.method, method, sizeof(gmx_ana_selmethod_t));
530     _gmx_selelem_init_method_params(sel, scanner);
531 }
532
533 /*! \brief
534  * Initializes the reference position calculation for a \ref SEL_EXPRESSION
535  * element.
536  *
537  * \param[in,out] pcc    Position calculation collection to use.
538  * \param[in,out] sel    Selection element to initialize.
539  * \param[in]     rpost  Reference position type to use (NULL = default).
540  * \param[in]     scanner Scanner data structure.
541  * \returns       0 on success, a non-zero error code on error.
542  */
543 static void
544 set_refpos_type(gmx::PositionCalculationCollection *pcc,
545                 const SelectionTreeElementPointer &sel,
546                 const char *rpost, yyscan_t scanner)
547 {
548     if (!rpost)
549     {
550         return;
551     }
552
553     if (sel->u.expr.method->pupdate)
554     {
555         /* By default, use whole residues/molecules. */
556         sel->u.expr.pc
557             = pcc->createCalculationFromEnum(rpost, POS_COMPLWHOLE);
558     }
559     else
560     {
561         // TODO: Should this be treated as a real error?
562         _gmx_selparser_error(scanner, "modifier '%s' is not applicable for '%s'",
563                              rpost, sel->u.expr.method->name);
564     }
565 }
566
567 gmx::SelectionTreeElementPointer
568 _gmx_sel_init_arithmetic(const gmx::SelectionTreeElementPointer &left,
569                          const gmx::SelectionTreeElementPointer &right,
570                          char op, yyscan_t /* scanner */)
571 {
572     SelectionTreeElementPointer sel(new SelectionTreeElement(SEL_ARITHMETIC));
573     sel->v.type        = REAL_VALUE;
574     switch (op)
575     {
576         case '+': sel->u.arith.type = ARITH_PLUS; break;
577         case '-': sel->u.arith.type = (right ? ARITH_MINUS : ARITH_NEG); break;
578         case '*': sel->u.arith.type = ARITH_MULT; break;
579         case '/': sel->u.arith.type = ARITH_DIV;  break;
580         case '^': sel->u.arith.type = ARITH_EXP;  break;
581     }
582     char               buf[2];
583     buf[0] = op;
584     buf[1] = 0;
585     sel->setName(buf);
586     sel->u.arith.opstr = strdup(buf);
587     sel->child         = left;
588     sel->child->next   = right;
589     return sel;
590 }
591
592 /*!
593  * \param[in]  left   Selection element for the left hand side.
594  * \param[in]  right  Selection element for the right hand side.
595  * \param[in]  cmpop  String representation of the comparison operator.
596  * \param[in]  scanner Scanner data structure.
597  * \returns    The created selection element.
598  *
599  * This function handles the creation of a gmx::SelectionTreeElement object for
600  * comparison expressions.
601  */
602 SelectionTreeElementPointer
603 _gmx_sel_init_comparison(const gmx::SelectionTreeElementPointer &left,
604                          const gmx::SelectionTreeElementPointer &right,
605                          const char *cmpop, yyscan_t scanner)
606 {
607     gmx::MessageStringCollector *errors = _gmx_sel_lexer_error_reporter(scanner);
608     gmx::MessageStringContext    context(errors, "In comparison initialization");
609
610     SelectionTreeElementPointer  sel(new SelectionTreeElement(SEL_EXPRESSION));
611     _gmx_selelem_set_method(sel, &sm_compare, scanner);
612
613     SelectionParserParameterList params;
614     const char                  *name;
615     // Create the parameter for the left expression.
616     name  = left->v.type == INT_VALUE ? "int1" : "real1";
617     params.push_back(SelectionParserParameter::createFromExpression(name, left));
618     // Create the parameter for the right expression.
619     name  = right->v.type == INT_VALUE ? "int2" : "real2";
620     params.push_back(SelectionParserParameter::createFromExpression(name, right));
621     // Create the parameter for the operator.
622     params.push_back(
623             SelectionParserParameter::create(
624                     "op", SelectionParserValue::createString(cmpop)));
625     if (!_gmx_sel_parse_params(params, sel->u.expr.method->nparams,
626                                sel->u.expr.method->param, sel, scanner))
627     {
628         return SelectionTreeElementPointer();
629     }
630
631     return sel;
632 }
633
634 /*! \brief
635  * Implementation method for keyword expression creation.
636  *
637  * \param[in]  method Method to use.
638  * \param[in]  matchType String matching type (only used if \p method is
639  *      a string keyword and \p args is not empty.
640  * \param[in]  args   Pointer to the first argument.
641  * \param[in]  rpost  Reference position type to use (NULL = default).
642  * \param[in]  scanner Scanner data structure.
643  * \returns    The created selection element.
644  *
645  * This function handles the creation of a gmx::SelectionTreeElement object for
646  * selection methods that do not take parameters.
647  */
648 static SelectionTreeElementPointer
649 init_keyword_internal(gmx_ana_selmethod_t *method,
650                       gmx::SelectionStringMatchType matchType,
651                       SelectionParserValueListPointer args,
652                       const char *rpost, yyscan_t scanner)
653 {
654     gmx_ana_selcollection_t     *sc = _gmx_sel_lexer_selcollection(scanner);
655
656     gmx::MessageStringCollector *errors = _gmx_sel_lexer_error_reporter(scanner);
657     char  buf[128];
658     sprintf(buf, "In keyword '%s'", method->name);
659     gmx::MessageStringContext  context(errors, buf);
660
661     if (method->nparams > 0)
662     {
663         // TODO: Would assert be better?
664         GMX_THROW(gmx::InternalError(
665                           "Keyword initialization called with non-keyword method"));
666     }
667
668     SelectionTreeElementPointer root(new SelectionTreeElement(SEL_EXPRESSION));
669     SelectionTreeElementPointer child = root;
670     _gmx_selelem_set_method(child, method, scanner);
671
672     /* Initialize the evaluation of keyword matching if values are provided */
673     if (args)
674     {
675         gmx_ana_selmethod_t *kwmethod;
676         switch (method->type)
677         {
678             case INT_VALUE:  kwmethod = &sm_keyword_int;  break;
679             case REAL_VALUE: kwmethod = &sm_keyword_real; break;
680             case STR_VALUE:  kwmethod = &sm_keyword_str;  break;
681             default:
682                 GMX_THROW(gmx::InternalError(
683                                   "Unknown type for keyword selection"));
684         }
685         /* Initialize the selection element */
686         root.reset(new SelectionTreeElement(SEL_EXPRESSION));
687         _gmx_selelem_set_method(root, kwmethod, scanner);
688         if (method->type == STR_VALUE)
689         {
690             _gmx_selelem_set_kwstr_match_type(root, matchType);
691         }
692         SelectionParserParameterList params;
693         params.push_back(
694                 SelectionParserParameter::createFromExpression(NULL, child));
695         params.push_back(SelectionParserParameter::create(NULL, move(args)));
696         if (!_gmx_sel_parse_params(params, root->u.expr.method->nparams,
697                                    root->u.expr.method->param, root, scanner))
698         {
699             return SelectionTreeElementPointer();
700         }
701     }
702     set_refpos_type(&sc->pcc, child, rpost, scanner);
703
704     return root;
705 }
706
707 /*!
708  * \param[in]  method Method to use.
709  * \param[in]  args   Pointer to the first argument.
710  * \param[in]  rpost  Reference position type to use (NULL = default).
711  * \param[in]  scanner Scanner data structure.
712  * \returns    The created selection element.
713  *
714  * This function handles the creation of a gmx::SelectionTreeElement object for
715  * selection methods that do not take parameters.
716  */
717 SelectionTreeElementPointer
718 _gmx_sel_init_keyword(gmx_ana_selmethod_t *method,
719                       gmx::SelectionParserValueListPointer args,
720                       const char *rpost, yyscan_t scanner)
721 {
722     return init_keyword_internal(method, gmx::eStringMatchType_Auto, move(args),
723                                  rpost, scanner);
724 }
725
726 /*!
727  * \param[in]  method    Method to use.
728  * \param[in]  matchType String matching type.
729  * \param[in]  args      Pointer to the first argument.
730  * \param[in]  rpost     Reference position type to use (NULL = default).
731  * \param[in]  scanner   Scanner data structure.
732  * \returns    The created selection element.
733  *
734  * This function handles the creation of a gmx::SelectionTreeElement object for
735  * keyword string matching.
736  */
737 SelectionTreeElementPointer
738 _gmx_sel_init_keyword_strmatch(gmx_ana_selmethod_t *method,
739                                gmx::SelectionStringMatchType matchType,
740                                gmx::SelectionParserValueListPointer args,
741                                const char *rpost, yyscan_t scanner)
742 {
743     GMX_RELEASE_ASSERT(method->type == STR_VALUE,
744                        "String keyword method called for a non-string-valued method");
745     GMX_RELEASE_ASSERT(args && !args->empty(),
746                        "String keyword matching method called without any values");
747     return init_keyword_internal(method, matchType, move(args), rpost, scanner);
748 }
749
750 /*!
751  * \param[in]  method Method to use for initialization.
752  * \param[in]  params Pointer to the first parameter.
753  * \param[in]  rpost  Reference position type to use (NULL = default).
754  * \param[in]  scanner Scanner data structure.
755  * \returns    The created selection element.
756  *
757  * This function handles the creation of a gmx::SelectionTreeElement object for
758  * selection methods that take parameters.
759  *
760  * Part of the behavior of the \c same selection keyword is hardcoded into
761  * this function (or rather, into _gmx_selelem_custom_init_same()) to allow the
762  * use of any keyword in \c "same KEYWORD as" without requiring special
763  * handling somewhere else (or sacrificing the simple syntax).
764  */
765 SelectionTreeElementPointer
766 _gmx_sel_init_method(gmx_ana_selmethod_t                      *method,
767                      gmx::SelectionParserParameterListPointer  params,
768                      const char *rpost, yyscan_t scanner)
769 {
770     gmx_ana_selcollection_t     *sc = _gmx_sel_lexer_selcollection(scanner);
771     int                          rc;
772
773     gmx::MessageStringCollector *errors = _gmx_sel_lexer_error_reporter(scanner);
774     char  buf[128];
775     sprintf(buf, "In keyword '%s'", method->name);
776     gmx::MessageStringContext  context(errors, buf);
777
778     _gmx_sel_finish_method(scanner);
779     /* The "same" keyword needs some custom massaging of the parameters. */
780     rc = _gmx_selelem_custom_init_same(&method, params, scanner);
781     if (rc != 0)
782     {
783         return SelectionTreeElementPointer();
784     }
785     SelectionTreeElementPointer root(new SelectionTreeElement(SEL_EXPRESSION));
786     _gmx_selelem_set_method(root, method, scanner);
787     /* Process the parameters */
788     if (!_gmx_sel_parse_params(*params, root->u.expr.method->nparams,
789                                root->u.expr.method->param, root, scanner))
790     {
791         return SelectionTreeElementPointer();
792     }
793     set_refpos_type(&sc->pcc, root, rpost, scanner);
794
795     return root;
796 }
797
798 /*!
799  * \param[in]  method Modifier to use for initialization.
800  * \param[in]  params Pointer to the first parameter.
801  * \param[in]  sel    Selection element that the modifier should act on.
802  * \param[in]  scanner Scanner data structure.
803  * \returns    The created selection element.
804  *
805  * This function handles the creation of a gmx::SelectionTreeElement object for
806  * selection modifiers.
807  */
808 SelectionTreeElementPointer
809 _gmx_sel_init_modifier(gmx_ana_selmethod_t                      *method,
810                        gmx::SelectionParserParameterListPointer  params,
811                        const gmx::SelectionTreeElementPointer   &sel,
812                        yyscan_t                                  scanner)
813 {
814     gmx::MessageStringCollector *errors = _gmx_sel_lexer_error_reporter(scanner);
815     char  buf[128];
816     sprintf(buf, "In keyword '%s'", method->name);
817     gmx::MessageStringContext  context(errors, buf);
818
819     _gmx_sel_finish_method(scanner);
820     SelectionTreeElementPointer modifier(new SelectionTreeElement(SEL_MODIFIER));
821     _gmx_selelem_set_method(modifier, method, scanner);
822     SelectionTreeElementPointer root;
823     if (method->type == NO_VALUE)
824     {
825         SelectionTreeElementPointer child = sel;
826         while (child->next)
827         {
828             child = child->next;
829         }
830         child->next = modifier;
831         root        = sel;
832     }
833     else
834     {
835         params->push_front(
836                 SelectionParserParameter::createFromExpression(NULL, sel));
837         root = modifier;
838     }
839     /* Process the parameters */
840     if (!_gmx_sel_parse_params(*params, modifier->u.expr.method->nparams,
841                                modifier->u.expr.method->param, modifier, scanner))
842     {
843         return SelectionTreeElementPointer();
844     }
845
846     return root;
847 }
848
849 /*!
850  * \param[in]  expr    Input selection element for the position calculation.
851  * \param[in]  type    Reference position type or NULL for default.
852  * \param[in]  scanner Scanner data structure.
853  * \returns    The created selection element.
854  *
855  * This function handles the creation of a gmx::SelectionTreeElement object for
856  * evaluation of reference positions.
857  */
858 SelectionTreeElementPointer
859 _gmx_sel_init_position(const gmx::SelectionTreeElementPointer &expr,
860                        const char *type, yyscan_t scanner)
861 {
862     gmx::MessageStringCollector *errors = _gmx_sel_lexer_error_reporter(scanner);
863     char  buf[128];
864     sprintf(buf, "In position evaluation");
865     gmx::MessageStringContext   context(errors, buf);
866
867     SelectionTreeElementPointer root(new SelectionTreeElement(SEL_EXPRESSION));
868     _gmx_selelem_set_method(root, &sm_keyword_pos, scanner);
869     _gmx_selelem_set_kwpos_type(root.get(), type);
870     /* Create the parameters for the parameter parser. */
871     SelectionParserParameterList params;
872     params.push_back(SelectionParserParameter::createFromExpression(NULL, expr));
873     /* Parse the parameters. */
874     if (!_gmx_sel_parse_params(params, root->u.expr.method->nparams,
875                                root->u.expr.method->param, root, scanner))
876     {
877         return SelectionTreeElementPointer();
878     }
879
880     return root;
881 }
882
883 /*!
884  * \param[in] x,y,z  Coordinates for the position.
885  * \returns   The creates selection element.
886  */
887 SelectionTreeElementPointer
888 _gmx_sel_init_const_position(real x, real y, real z)
889 {
890     rvec                        pos;
891
892     SelectionTreeElementPointer sel(new SelectionTreeElement(SEL_CONST));
893     _gmx_selelem_set_vtype(sel, POS_VALUE);
894     _gmx_selvalue_reserve(&sel->v, 1);
895     pos[XX] = x;
896     pos[YY] = y;
897     pos[ZZ] = z;
898     gmx_ana_pos_init_const(sel->v.u.p, pos);
899     return sel;
900 }
901
902 /*!
903  * \param[in] name  Name of an index group to search for.
904  * \param[in] scanner Scanner data structure.
905  * \returns   The created selection element.
906  *
907  * See gmx_ana_indexgrps_find() for information on how \p name is matched
908  * against the index groups.
909  */
910 SelectionTreeElementPointer
911 _gmx_sel_init_group_by_name(const char *name, yyscan_t scanner)
912 {
913
914     SelectionTreeElementPointer sel(new SelectionTreeElement(SEL_GROUPREF));
915     _gmx_selelem_set_vtype(sel, GROUP_VALUE);
916     sel->setName(gmx::formatString("group \"%s\"", name));
917     sel->u.gref.name = strdup(name);
918     sel->u.gref.id   = -1;
919
920     if (_gmx_sel_lexer_has_groups_set(scanner))
921     {
922         gmx_ana_indexgrps_t *grps = _gmx_sel_lexer_indexgrps(scanner);
923         sel->resolveIndexGroupReference(grps);
924     }
925
926     return sel;
927 }
928
929 /*!
930  * \param[in] id    Zero-based index number of the group to extract.
931  * \param[in] scanner Scanner data structure.
932  * \returns   The created selection element.
933  */
934 SelectionTreeElementPointer
935 _gmx_sel_init_group_by_id(int id, yyscan_t scanner)
936 {
937     SelectionTreeElementPointer sel(new SelectionTreeElement(SEL_GROUPREF));
938     _gmx_selelem_set_vtype(sel, GROUP_VALUE);
939     sel->setName(gmx::formatString("group %d", id));
940     sel->u.gref.name = NULL;
941     sel->u.gref.id   = id;
942
943     if (_gmx_sel_lexer_has_groups_set(scanner))
944     {
945         gmx_ana_indexgrps_t *grps = _gmx_sel_lexer_indexgrps(scanner);
946         sel->resolveIndexGroupReference(grps);
947     }
948
949     return sel;
950 }
951
952 /*!
953  * \param[in,out] sel  Value of the variable.
954  * \returns       The created selection element that references \p sel.
955  *
956  * The reference count of \p sel is updated, but no other modifications are
957  * made.
958  */
959 SelectionTreeElementPointer
960 _gmx_sel_init_variable_ref(const gmx::SelectionTreeElementPointer &sel)
961 {
962     SelectionTreeElementPointer ref;
963
964     if (sel->v.type == POS_VALUE && sel->type == SEL_CONST)
965     {
966         ref = sel;
967     }
968     else
969     {
970         ref.reset(new SelectionTreeElement(SEL_SUBEXPRREF));
971         _gmx_selelem_set_vtype(ref, sel->v.type);
972         ref->setName(sel->name());
973         ref->child = sel;
974     }
975     return ref;
976 }
977
978 /*!
979  * \param[in]  name     Name for the selection
980  *     (if NULL, a default name is constructed).
981  * \param[in]  sel      The selection element that evaluates the selection.
982  * \param      scanner  Scanner data structure.
983  * \returns    The created root selection element.
984  *
985  * This function handles the creation of root (\ref SEL_ROOT)
986  * gmx::SelectionTreeElement objects for selections.
987  */
988 SelectionTreeElementPointer
989 _gmx_sel_init_selection(const char                             *name,
990                         const gmx::SelectionTreeElementPointer &sel,
991                         yyscan_t                                scanner)
992 {
993     gmx::MessageStringCollector *errors = _gmx_sel_lexer_error_reporter(scanner);
994     char  buf[1024];
995     sprintf(buf, "In selection '%s'", _gmx_sel_lexer_pselstr(scanner));
996     gmx::MessageStringContext  context(errors, buf);
997
998     if (sel->v.type != POS_VALUE)
999     {
1000         /* FIXME: Better handling of this error */
1001         GMX_THROW(gmx::InternalError(
1002                           "Each selection must evaluate to a position"));
1003     }
1004
1005     SelectionTreeElementPointer root(new SelectionTreeElement(SEL_ROOT));
1006     root->child = sel;
1007     if (name)
1008     {
1009         root->setName(name);
1010     }
1011     /* Update the flags */
1012     _gmx_selelem_update_flags(root, scanner);
1013
1014     root->fillNameIfMissing(_gmx_sel_lexer_pselstr(scanner));
1015
1016     /* Print out some information if the parser is interactive */
1017     if (_gmx_sel_is_lexer_interactive(scanner))
1018     {
1019         fprintf(stderr, "Selection '%s' parsed\n",
1020                 _gmx_sel_lexer_pselstr(scanner));
1021     }
1022
1023     return root;
1024 }
1025
1026
1027 /*!
1028  * \param[in]  name     Name of the variable.
1029  * \param[in]  expr     The selection element that evaluates the variable.
1030  * \param      scanner  Scanner data structure.
1031  * \returns    The created root selection element.
1032  *
1033  * This function handles the creation of root gmx::SelectionTreeElement objects
1034  * for variable assignments. A \ref SEL_ROOT element and a \ref SEL_SUBEXPR
1035  * element are both created.
1036  */
1037 SelectionTreeElementPointer
1038 _gmx_sel_assign_variable(const char                             *name,
1039                          const gmx::SelectionTreeElementPointer &expr,
1040                          yyscan_t                                scanner)
1041 {
1042     gmx_ana_selcollection_t     *sc      = _gmx_sel_lexer_selcollection(scanner);
1043     const char                  *pselstr = _gmx_sel_lexer_pselstr(scanner);
1044     SelectionTreeElementPointer  root;
1045
1046     gmx::MessageStringCollector *errors = _gmx_sel_lexer_error_reporter(scanner);
1047     char  buf[1024];
1048     sprintf(buf, "In selection '%s'", pselstr);
1049     gmx::MessageStringContext  context(errors, buf);
1050
1051     _gmx_selelem_update_flags(expr, scanner);
1052     /* Check if this is a constant non-group value */
1053     if (expr->type == SEL_CONST && expr->v.type != GROUP_VALUE)
1054     {
1055         /* If so, just assign the constant value to the variable */
1056         sc->symtab->addVariable(name, expr);
1057         goto finish;
1058     }
1059     /* Check if we are assigning a variable to another variable */
1060     if (expr->type == SEL_SUBEXPRREF)
1061     {
1062         /* If so, make a simple alias */
1063         sc->symtab->addVariable(name, expr->child);
1064         goto finish;
1065     }
1066     /* Create the root element */
1067     root.reset(new SelectionTreeElement(SEL_ROOT));
1068     root->setName(name);
1069     /* Create the subexpression element */
1070     root->child.reset(new SelectionTreeElement(SEL_SUBEXPR));
1071     root->child->setName(name);
1072     _gmx_selelem_set_vtype(root->child, expr->v.type);
1073     root->child->child  = expr;
1074     /* Update flags */
1075     _gmx_selelem_update_flags(root, scanner);
1076     /* Add the variable to the symbol table */
1077     sc->symtab->addVariable(name, root->child);
1078 finish:
1079     srenew(sc->varstrs, sc->nvars + 1);
1080     sc->varstrs[sc->nvars] = strdup(pselstr);
1081     ++sc->nvars;
1082     if (_gmx_sel_is_lexer_interactive(scanner))
1083     {
1084         fprintf(stderr, "Variable '%s' parsed\n", pselstr);
1085     }
1086     return root;
1087 }
1088
1089 /*!
1090  * \param         sel   Selection to append (can be NULL, in which
1091  *   case nothing is done).
1092  * \param         last  Last selection, or NULL if not present or not known.
1093  * \param         scanner  Scanner data structure.
1094  * \returns       The last selection after the append.
1095  *
1096  * Appends \p sel after the last root element, and returns either \p sel
1097  * (if it was non-NULL) or the last element (if \p sel was NULL).
1098  */
1099 SelectionTreeElementPointer
1100 _gmx_sel_append_selection(const gmx::SelectionTreeElementPointer &sel,
1101                           gmx::SelectionTreeElementPointer        last,
1102                           yyscan_t                                scanner)
1103 {
1104     gmx_ana_selcollection_t *sc = _gmx_sel_lexer_selcollection(scanner);
1105
1106     /* Append sel after last, or the last element of sc if last is NULL */
1107     if (last)
1108     {
1109         last->next = sel;
1110     }
1111     else
1112     {
1113         if (sc->root)
1114         {
1115             last = sc->root;
1116             while (last->next)
1117             {
1118                 last = last->next;
1119             }
1120             last->next = sel;
1121         }
1122         else
1123         {
1124             sc->root = sel;
1125         }
1126     }
1127     /* Initialize a selection object if necessary */
1128     if (sel)
1129     {
1130         last = sel;
1131         /* Add the new selection to the collection if it is not a variable. */
1132         if (sel->child->type != SEL_SUBEXPR)
1133         {
1134             gmx::SelectionDataPointer selPtr(
1135                     new gmx::internal::SelectionData(
1136                             sel.get(), _gmx_sel_lexer_pselstr(scanner)));
1137             sc->sel.push_back(gmx::move(selPtr));
1138         }
1139     }
1140     /* Clear the selection string now that we've saved it */
1141     _gmx_sel_lexer_clear_pselstr(scanner);
1142     return last;
1143 }
1144
1145 /*!
1146  * \param[in] scanner Scanner data structure.
1147  * \returns   true if the parser should finish, false if parsing should
1148  *   continue.
1149  *
1150  * This function is called always after _gmx_sel_append_selection() to
1151  * check whether a sufficient number of selections has already been provided.
1152  * This is used to terminate interactive parsers when the correct number of
1153  * selections has been provided.
1154  */
1155 bool
1156 _gmx_sel_parser_should_finish(yyscan_t scanner)
1157 {
1158     gmx_ana_selcollection_t *sc = _gmx_sel_lexer_selcollection(scanner);
1159     return (int)sc->sel.size() == _gmx_sel_lexer_exp_selcount(scanner);
1160 }