05ea9b6c6f1f65b4d9f3eeb06fbabb9801a8fcb9
[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,2015,2016,2017,2018,2019, 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::parseInteractive(),
70  *    SelectionCollection::parseFromStdin(),
71  *    SelectionCollection::parseFromFile(), and
72  *    SelectionCollection::parseFromString().
73  *
74  * The basic control flow in the parser is as follows: when a parser function
75  * in SelectionCollection gets called, it performs some
76  * initialization, and then calls the _gmx_sel_yyparse() function generated
77  * by Bison. This function then calls _gmx_sel_yylex() to repeatedly read
78  * tokens from the input (more complex tasks related to token recognition
79  * and bookkeeping are done by functions in scanner_internal.cpp) and uses the
80  * grammar rules to decide what to do with them. Whenever a grammar rule
81  * matches, a corresponding function in parsetree.cpp is called to construct
82  * either a temporary representation for the object or a
83  * gmx::SelectionTreeElement object
84  * (some simple rules are handled internally in parser.y).
85  * When a complete selection has been parsed, the functions in parsetree.cpp
86  * also take care of updating the ::gmx_ana_selcollection_t structure
87  * appropriately.
88  *
89  * The rest of this page describes the resulting gmx::SelectionTreeElement
90  * object tree.
91  * Before the selections can be evaluated, this tree needs to be passed to
92  * the selection compiler, which is described on a separate page:
93  * \ref page_module_selection_compiler
94  *
95  *
96  * \section selparser_tree Element tree constructed by the parser
97  *
98  * The parser initializes the following fields in all selection elements:
99  * gmx::SelectionTreeElement::name, gmx::SelectionTreeElement::type,
100  * gmx::SelectionTreeElement::v\c .type,
101  * gmx::SelectionTreeElement::flags, gmx::SelectionTreeElement::child, and
102  * gmx::SelectionTreeElement::next.
103  * Some other fields are also initialized for particular element types as
104  * discussed below.
105  * Fields that are not initialized are set to zero, NULL, or other similar
106  * value.
107  *
108  *
109  * \subsection selparser_tree_root Root elements
110  *
111  * The parser creates a \ref SEL_ROOT selection element for each variable
112  * assignment and each selection. However, there are two exceptions that do
113  * not result in a \ref SEL_ROOT element (in these cases, only the symbol
114  * table is modified):
115  *  - Variable assignments that assign a variable to another variable.
116  *  - Variable assignments that assign a non-group constant.
117  *  .
118  * The \ref SEL_ROOT elements are linked together in a chain in the same order
119  * as in the input.
120  *
121  * The children of the \ref SEL_ROOT elements can be used to distinguish
122  * the two types of root elements from each other:
123  *  - For variable assignments, the first and only child is always
124  *    a \ref SEL_SUBEXPR element.
125  *  - For selections, the first child is a \ref SEL_EXPRESSION or a
126  *    \ref SEL_MODIFIER element that evaluates the final positions (if the
127  *    selection defines a constant position, the child is a \ref SEL_CONST).
128  *    The rest of the children are \ref SEL_MODIFIER elements with
129  *    \ref NO_VALUE, in the order given by the user.
130  *  .
131  * The name of the selection/variable is stored in
132  * gmx::SelectionTreeElement::cgrp\c .name.
133  * It is set to either the name provided by the user or the selection string
134  * for selections not explicitly named by the user.
135  * \ref SEL_ROOT or \ref SEL_SUBEXPR elements do not appear anywhere else.
136  *
137  *
138  * \subsection selparser_tree_const Constant elements
139  *
140  * \ref SEL_CONST elements are created for every constant that is required
141  * for later evaluation.
142  * Currently, \ref SEL_CONST elements can be present for
143  *  - selections that consist of a constant position,
144  *  - \ref GROUP_VALUE method parameters if provided using external index
145  *    groups,
146  *  .
147  * For group-valued elements, the value is stored in
148  * gmx::SelectionTreeElement::cgrp; other types of values are stored in
149  * gmx::SelectionTreeElement::v.
150  * Constants that appear as parameters for selection methods are not present
151  * in the selection tree unless they have \ref GROUP_VALUE.
152  * \ref SEL_CONST elements have no children.
153  *
154  *
155  * \subsection selparser_tree_method Method evaluation elements
156  *
157  * \ref SEL_EXPRESSION and \ref SEL_MODIFIER elements are treated very
158  * similarly. The \c gmx_ana_selmethod_t structure corresponding to the
159  * evaluation method is in gmx::SelectionTreeElement::method, and the method
160  * data in gmx::SelectionTreeElement::mdata has been allocated using
161  * sel_datafunc().
162  * If a non-standard reference position type was set,
163  * gmx::SelectionTreeElement::pc has also been created, but only the type has
164  * been set.
165  * All children of these elements are of the type \ref SEL_SUBEXPRREF, and
166  * each describes a selection that needs to be evaluated to obtain a value
167  * for one parameter of the method.
168  * No children are present for parameters that were given a constant
169  * non-\ref GROUP_VALUE value.
170  * The children are sorted in the order in which the parameters appear in the
171  * \ref gmx_ana_selmethod_t structure.
172  *
173  * In addition to actual selection keywords, \ref SEL_EXPRESSION elements
174  * are used internally to implement numerical comparisons (e.g., "x < 5")
175  * and keyword matching (e.g., "resnr 1 to 3" or "name CA").
176  *
177  *
178  * \subsection selparser_tree_subexpr Subexpression elements
179  *
180  * \ref SEL_SUBEXPR elements only appear for variables, as described above.
181  * gmx::SelectionTreeElement::name points to the name of the variable (from the
182  * \ref SEL_ROOT element).
183  * The element always has exactly one child, which represents the value of
184  * the variable.
185  *
186  * \ref SEL_SUBEXPRREF elements are used for two purposes:
187  *  - Variable references that need to be evaluated (i.e., there is a
188  *    \ref SEL_SUBEXPR element for the variable) are represented using
189  *    \ref SEL_SUBEXPRREF elements.
190  *    In this case, gmx::SelectionTreeElement::param is NULL, and the first and
191  *    only child of the element is the \ref SEL_SUBEXPR element of the
192  *    variable.
193  *    Such references can appear anywhere where the variable value
194  *    (the child of the \ref SEL_SUBEXPR element) would be valid.
195  *  - Children of \ref SEL_EXPRESSION and \ref SEL_MODIFIER elements are
196  *    always of this type. For these elements, gmx::SelectionTreeElement::param
197  *    is initialized to point to the parameter that receives the value from
198  *    the expression.
199  *    Each such element has exactly one child, which can be of any type;
200  *    the \ref SEL_SUBEXPR element of a variable is used if the value comes
201  *    from a variable, otherwise the child type is not \ref SEL_SUBEXPR.
202  *
203  *
204  * \subsection selparser_tree_bool Boolean elements
205  *
206  * One \ref SEL_BOOLEAN element is created for each boolean keyword in the
207  * input, and the tree structure represents the evaluation order.
208  * The gmx::SelectionTreeElement::boolt type gives the type of the operation.
209  * Each element has exactly two children (one for \ref BOOL_NOT elements),
210  * which are in the order given in the input.
211  * The children always have \ref GROUP_VALUE, but different element types
212  * are possible.
213  *
214  *
215  * \subsection selparser_tree_arith Arithmetic elements
216  *
217  * One \ref SEL_ARITHMETIC element is created for each arithmetic operation in
218  * the input, and the tree structure represents the evaluation order.
219  * The gmx::SelectionTreeElement::optype type gives the name of the operation.
220  * Each element has exactly two children (one for unary negation elements),
221  * which are in the order given in the input.
222  */
223 #include "gmxpre.h"
224
225 #include "parsetree.h"
226
227 #include <cstdarg>
228 #include <cstdio>
229
230 #include <exception>
231 #include <memory>
232
233 #include "gromacs/selection/selection.h"
234 #include "gromacs/utility/cstringutil.h"
235 #include "gromacs/utility/exceptions.h"
236 #include "gromacs/utility/smalloc.h"
237 #include "gromacs/utility/stringutil.h"
238 #include "gromacs/utility/textwriter.h"
239
240 #include "keywords.h"
241 #include "poscalc.h"
242 #include "scanner.h"
243 #include "selectioncollection_impl.h"
244 #include "selelem.h"
245 #include "selmethod.h"
246 #include "symrec.h"
247
248 using gmx::SelectionLocation;
249 using gmx::SelectionParserParameter;
250 using gmx::SelectionParserParameterList;
251 using gmx::SelectionParserParameterListPointer;
252 using gmx::SelectionParserValue;
253 using gmx::SelectionParserValueList;
254 using gmx::SelectionParserValueListPointer;
255 using gmx::SelectionTreeElement;
256 using gmx::SelectionTreeElementPointer;
257
258 namespace
259 {
260
261 /*! \brief
262  * Formats context string for errors.
263  *
264  * The returned string is used as the context for errors reported during
265  * parsing.
266  */
267 std::string formatCurrentErrorContext(yyscan_t scanner)
268 {
269     return gmx::formatString("While parsing '%s'", _gmx_sel_lexer_get_current_text(scanner).c_str());
270 }
271
272 } // namespace
273
274 bool _gmx_selparser_handle_exception(yyscan_t scanner, std::exception* ex)
275 {
276     try
277     {
278         bool                   canContinue      = false;
279         gmx::GromacsException* gromacsException = dynamic_cast<gmx::GromacsException*>(ex);
280         if (gromacsException != nullptr)
281         {
282             gromacsException->prependContext(formatCurrentErrorContext(scanner));
283             canContinue = (dynamic_cast<gmx::UserInputError*>(ex) != nullptr);
284         }
285         _gmx_sel_lexer_set_exception(scanner, std::current_exception());
286         return canContinue;
287     }
288     catch (const std::exception&)
289     {
290         _gmx_sel_lexer_set_exception(scanner, std::current_exception());
291         return false;
292     }
293 }
294
295 bool _gmx_selparser_handle_error(yyscan_t scanner)
296 {
297     std::string context(gmx::formatString("Invalid selection '%s'", _gmx_sel_lexer_pselstr(scanner)));
298     // The only way to prepend context to the exception is to rethrow it.
299     try
300     {
301         _gmx_sel_lexer_rethrow_exception_if_occurred(scanner);
302     }
303     catch (gmx::UserInputError& ex)
304     {
305         ex.prependContext(context);
306         gmx::TextWriter* statusWriter = _gmx_sel_lexer_get_status_writer(scanner);
307         if (statusWriter != nullptr)
308         {
309             gmx::formatExceptionMessageToWriter(statusWriter, ex);
310             return true;
311         }
312         throw;
313     }
314     catch (gmx::GromacsException& ex)
315     {
316         ex.prependContext(context);
317         throw;
318     }
319     GMX_RELEASE_ASSERT(false, "All parsing errors should result in a captured exception");
320     return false; // Some compilers will not believe that the above never returns.
321 }
322
323 namespace gmx
324 {
325
326 /********************************************************************
327  * SelectionParserValue
328  */
329
330 SelectionParserValue::SelectionParserValue(e_selvalue_t type, const SelectionLocation& location) :
331     type(type),
332     location_(location)
333 {
334     memset(&u, 0, sizeof(u));
335 }
336
337 SelectionParserValue::SelectionParserValue(const SelectionTreeElementPointer& expr) :
338     type(expr->v.type),
339     expr(expr),
340     location_(expr->location())
341 {
342     memset(&u, 0, sizeof(u));
343 }
344
345 /********************************************************************
346  * SelectionParserParameter
347  */
348
349 SelectionParserParameter::SelectionParserParameter(const char*                     name,
350                                                    SelectionParserValueListPointer values,
351                                                    const SelectionLocation&        location) :
352     name_(name != nullptr ? name : ""),
353     location_(location),
354     values_(values ? std::move(values) : std::make_unique<SelectionParserValueList>())
355 {
356 }
357
358 } // namespace gmx
359
360 /*!
361  * \param[in,out] sel  Root of the selection element tree to initialize.
362  *
363  * Propagates the \ref SEL_DYNAMIC flag from the children of \p sel to \p sel
364  * (if any child of \p sel is dynamic, \p sel is also marked as such).
365  * The \ref SEL_DYNAMIC flag is also set for \ref SEL_EXPRESSION elements with
366  * a dynamic method.
367  * Also, sets one of the \ref SEL_SINGLEVAL, \ref SEL_ATOMVAL, or
368  * \ref SEL_VARNUMVAL flags, either based on the children or on the type of
369  * the selection method.
370  * If the types of the children conflict, an error is returned.
371  *
372  * The flags of the children of \p sel are also updated if not done earlier.
373  * The flags are initialized only once for any element; if \ref SEL_FLAGSSET
374  * is set for an element, the function returns immediately, and the recursive
375  * operation does not descend beyond such elements.
376  */
377 void _gmx_selelem_update_flags(const gmx::SelectionTreeElementPointer& sel)
378 {
379     bool bUseChildType = false;
380     bool bOnlySingleChildren;
381
382     /* Return if the flags have already been set */
383     if (sel->flags & SEL_FLAGSSET)
384     {
385         return;
386     }
387     /* Set the flags based on the current element type */
388     switch (sel->type)
389     {
390         case SEL_CONST:
391         case SEL_GROUPREF:
392             sel->flags |= SEL_SINGLEVAL;
393             bUseChildType = false;
394             break;
395
396         case SEL_EXPRESSION:
397             if (sel->u.expr.method->flags & SMETH_DYNAMIC)
398             {
399                 sel->flags |= SEL_DYNAMIC;
400             }
401             if (sel->u.expr.method->flags & SMETH_SINGLEVAL)
402             {
403                 sel->flags |= SEL_SINGLEVAL;
404             }
405             else if (sel->u.expr.method->flags & SMETH_VARNUMVAL)
406             {
407                 sel->flags |= SEL_VARNUMVAL;
408             }
409             else
410             {
411                 sel->flags |= SEL_ATOMVAL;
412             }
413             bUseChildType = false;
414             break;
415
416         case SEL_ARITHMETIC:
417             sel->flags |= SEL_ATOMVAL;
418             bUseChildType = false;
419             break;
420
421         case SEL_MODIFIER:
422             if (sel->v.type != NO_VALUE)
423             {
424                 sel->flags |= SEL_VARNUMVAL;
425             }
426             bUseChildType = false;
427             break;
428
429         case SEL_ROOT: bUseChildType = false; break;
430
431         case SEL_BOOLEAN:
432         case SEL_SUBEXPR:
433         case SEL_SUBEXPRREF: bUseChildType = true; break;
434     }
435     /* Loop through children to propagate their flags upwards */
436     bOnlySingleChildren               = true;
437     SelectionTreeElementPointer child = sel->child;
438     while (child)
439     {
440         /* Update the child */
441         _gmx_selelem_update_flags(child);
442         /* Propagate the dynamic and unsorted flags */
443         sel->flags |= (child->flags & (SEL_DYNAMIC | SEL_UNSORTED));
444         /* Propagate the type flag if necessary and check for problems */
445         if (bUseChildType)
446         {
447             if ((sel->flags & SEL_VALTYPEMASK) && !(sel->flags & child->flags & SEL_VALTYPEMASK))
448             {
449                 // TODO: Recollect when this is triggered, and whether the type
450                 // is appropriate.
451                 GMX_THROW(gmx::InvalidInputError("Invalid combination of selection expressions"));
452             }
453             sel->flags |= (child->flags & SEL_VALTYPEMASK);
454         }
455         if (!(child->flags & SEL_SINGLEVAL))
456         {
457             bOnlySingleChildren = false;
458         }
459
460         child = child->next;
461     }
462     /* For arithmetic expressions consisting only of single values,
463      * the result is also a single value. */
464     if (sel->type == SEL_ARITHMETIC && bOnlySingleChildren)
465     {
466         sel->flags = (sel->flags & ~SEL_VALTYPEMASK) | SEL_SINGLEVAL;
467     }
468     /* For root elements, the type should be propagated here, after the
469      * children have been updated. */
470     if (sel->type == SEL_ROOT)
471     {
472         GMX_ASSERT(sel->child, "Root elements should always have a child");
473         sel->flags |= (sel->child->flags & SEL_VALTYPEMASK);
474     }
475     /* Mark that the flags are set */
476     sel->flags |= SEL_FLAGSSET;
477 }
478
479 /*!
480  * \param[in,out] sel    Selection element to initialize.
481  * \param[in]     scanner Scanner data structure.
482  *
483  * A deep copy of the parameters is made to allow several
484  * expressions with the same method to coexist peacefully.
485  * Calls sel_datafunc() if one is specified for the method.
486  */
487 void _gmx_selelem_init_method_params(const gmx::SelectionTreeElementPointer& sel, yyscan_t scanner)
488 {
489     int                 nparams;
490     gmx_ana_selparam_t* orgparam;
491     gmx_ana_selparam_t* param;
492     int                 i;
493     void*               mdata;
494
495     nparams  = sel->u.expr.method->nparams;
496     orgparam = sel->u.expr.method->param;
497     snew(param, nparams);
498     memcpy(param, orgparam, nparams * sizeof(gmx_ana_selparam_t));
499     for (i = 0; i < nparams; ++i)
500     {
501         param[i].flags &= ~SPAR_SET;
502         _gmx_selvalue_setstore(&param[i].val, nullptr);
503         if (param[i].flags & SPAR_VARNUM)
504         {
505             param[i].val.nr = -1;
506         }
507         /* Duplicate the enum value array if it is given statically */
508         if ((param[i].flags & SPAR_ENUMVAL) && orgparam[i].val.u.ptr != nullptr)
509         {
510             int n;
511
512             /* Count the values */
513             n = 1;
514             while (orgparam[i].val.u.s[n] != nullptr)
515             {
516                 ++n;
517             }
518             _gmx_selvalue_reserve(&param[i].val, n + 1);
519             memcpy(param[i].val.u.s, orgparam[i].val.u.s, (n + 1) * sizeof(param[i].val.u.s[0]));
520         }
521     }
522     mdata = nullptr;
523     if (sel->u.expr.method->init_data)
524     {
525         mdata = sel->u.expr.method->init_data(nparams, param);
526     }
527     if (sel->u.expr.method->set_poscoll)
528     {
529         gmx_ana_selcollection_t* sc = _gmx_sel_lexer_selcollection(scanner);
530
531         sel->u.expr.method->set_poscoll(&sc->pcc, mdata);
532     }
533     /* Store the values */
534     sel->u.expr.method->param = param;
535     sel->u.expr.mdata         = mdata;
536 }
537
538 /*!
539  * \param[in,out] sel    Selection element to initialize.
540  * \param[in]     method Selection method to set.
541  * \param[in]     scanner Scanner data structure.
542  *
543  * Makes a copy of \p method and stores it in \p sel->u.expr.method,
544  * and calls _gmx_selelem_init_method_params();
545  */
546 void _gmx_selelem_set_method(const gmx::SelectionTreeElementPointer& sel,
547                              gmx_ana_selmethod_t*                    method,
548                              yyscan_t                                scanner)
549 {
550     _gmx_selelem_set_vtype(sel, method->type);
551     sel->setName(method->name);
552     snew(sel->u.expr.method, 1);
553     memcpy(sel->u.expr.method, method, sizeof(gmx_ana_selmethod_t));
554     _gmx_selelem_init_method_params(sel, scanner);
555 }
556
557 /*! \brief
558  * Initializes the reference position calculation for a \ref SEL_EXPRESSION
559  * element.
560  *
561  * \param[in,out] pcc    Position calculation collection to use.
562  * \param[in,out] sel    Selection element to initialize.
563  * \param[in]     rpost  Reference position type to use (NULL = default).
564  */
565 static void set_refpos_type(gmx::PositionCalculationCollection* pcc,
566                             const SelectionTreeElementPointer&  sel,
567                             const char*                         rpost)
568 {
569     if (!rpost)
570     {
571         return;
572     }
573
574     if (sel->u.expr.method->pupdate)
575     {
576         /* By default, use whole residues/molecules. */
577         sel->u.expr.pc = pcc->createCalculationFromEnum(rpost, POS_COMPLWHOLE);
578     }
579     else
580     {
581         std::string message = gmx::formatString(
582                 "Position modifiers ('%s') is not applicable for '%s'", rpost, sel->u.expr.method->name);
583         GMX_THROW(gmx::InvalidInputError(message));
584     }
585 }
586
587 gmx::SelectionTreeElementPointer _gmx_sel_init_arithmetic(const gmx::SelectionTreeElementPointer& left,
588                                                           const gmx::SelectionTreeElementPointer& right,
589                                                           char     op,
590                                                           yyscan_t scanner)
591 {
592     SelectionTreeElementPointer sel(
593             new SelectionTreeElement(SEL_ARITHMETIC, _gmx_sel_lexer_get_current_location(scanner)));
594     sel->v.type = REAL_VALUE;
595     switch (op)
596     {
597         case '+': sel->u.arith.type = ARITH_PLUS; break;
598         case '-': sel->u.arith.type = (right ? ARITH_MINUS : ARITH_NEG); break;
599         case '*': sel->u.arith.type = ARITH_MULT; break;
600         case '/': sel->u.arith.type = ARITH_DIV; break;
601         case '^': sel->u.arith.type = ARITH_EXP; break;
602     }
603     char buf[2]{ op, 0 };
604     sel->setName(buf);
605     sel->u.arith.opstr = gmx_strdup(buf);
606     sel->child         = left;
607     sel->child->next   = right;
608     return sel;
609 }
610
611 /*!
612  * \param[in]  left   Selection element for the left hand side.
613  * \param[in]  right  Selection element for the right hand side.
614  * \param[in]  cmpop  String representation of the comparison operator.
615  * \param[in]  scanner Scanner data structure.
616  * \returns    The created selection element.
617  *
618  * This function handles the creation of a gmx::SelectionTreeElement object for
619  * comparison expressions.
620  */
621 SelectionTreeElementPointer _gmx_sel_init_comparison(const gmx::SelectionTreeElementPointer& left,
622                                                      const gmx::SelectionTreeElementPointer& right,
623                                                      const char*                             cmpop,
624                                                      yyscan_t scanner)
625 {
626     SelectionTreeElementPointer sel(
627             new SelectionTreeElement(SEL_EXPRESSION, _gmx_sel_lexer_get_current_location(scanner)));
628     _gmx_selelem_set_method(sel, &sm_compare, scanner);
629
630     SelectionParserParameterList params;
631     const char*                  name;
632     // Create the parameter for the left expression.
633     name = left->v.type == INT_VALUE ? "int1" : "real1";
634     params.push_back(SelectionParserParameter::createFromExpression(name, left));
635     // Create the parameter for the right expression.
636     name = right->v.type == INT_VALUE ? "int2" : "real2";
637     params.push_back(SelectionParserParameter::createFromExpression(name, right));
638     // Create the parameter for the operator.
639     // TODO: Consider whether a proper location is needed.
640     SelectionLocation location(SelectionLocation::createEmpty());
641     params.push_back(SelectionParserParameter::create(
642             "op", SelectionParserValue::createString(cmpop, location), location));
643     _gmx_sel_parse_params(params, sel->u.expr.method->nparams, sel->u.expr.method->param, sel, scanner);
644
645     return sel;
646 }
647
648 /*! \brief
649  * Implementation method for keyword expression creation.
650  *
651  * \param[in]  method Method to use.
652  * \param[in]  matchType String matching type (only used if \p method is
653  *      a string keyword and \p args is not empty.
654  * \param[in]  args   Pointer to the first argument.
655  * \param[in]  rpost  Reference position type to use (NULL = default).
656  * \param[in]  scanner Scanner data structure.
657  * \returns    The created selection element.
658  *
659  * This function handles the creation of a gmx::SelectionTreeElement object for
660  * selection methods that do not take parameters.
661  */
662 static SelectionTreeElementPointer init_keyword_internal(gmx_ana_selmethod_t*            method,
663                                                          gmx::SelectionStringMatchType   matchType,
664                                                          SelectionParserValueListPointer args,
665                                                          const char*                     rpost,
666                                                          yyscan_t                        scanner)
667 {
668     gmx_ana_selcollection_t* sc = _gmx_sel_lexer_selcollection(scanner);
669
670     if (method->nparams > 0)
671     {
672         // TODO: Would assert be better?
673         GMX_THROW(gmx::InternalError("Keyword initialization called with non-keyword method"));
674     }
675
676     const SelectionLocation& location = _gmx_sel_lexer_get_current_location(scanner);
677     // TODO: If there are arguments, the location would be better as just the
678     // location of the keyword itself.
679     SelectionTreeElementPointer root(new SelectionTreeElement(SEL_EXPRESSION, location));
680     SelectionTreeElementPointer child = root;
681     _gmx_selelem_set_method(child, method, scanner);
682
683     /* Initialize the evaluation of keyword matching if values are provided */
684     if (args)
685     {
686         gmx_ana_selmethod_t* kwmethod;
687         switch (method->type)
688         {
689             case INT_VALUE: kwmethod = &sm_keyword_int; break;
690             case REAL_VALUE: kwmethod = &sm_keyword_real; break;
691             case STR_VALUE: kwmethod = &sm_keyword_str; break;
692             default: GMX_THROW(gmx::InternalError("Unknown type for keyword selection"));
693         }
694         /* Initialize the selection element */
695         root = std::make_shared<SelectionTreeElement>(SEL_EXPRESSION, location);
696         _gmx_selelem_set_method(root, kwmethod, scanner);
697         if (method->type == STR_VALUE)
698         {
699             _gmx_selelem_set_kwstr_match_type(root, matchType);
700         }
701         SelectionParserParameterList params;
702         params.push_back(SelectionParserParameter::createFromExpression(nullptr, child));
703         params.push_back(SelectionParserParameter::create(nullptr, std::move(args), location));
704         _gmx_sel_parse_params(params, root->u.expr.method->nparams, root->u.expr.method->param,
705                               root, scanner);
706     }
707     set_refpos_type(&sc->pcc, child, rpost);
708
709     return root;
710 }
711
712 /*!
713  * \param[in]  method Method to use.
714  * \param[in]  args   Pointer to the first argument.
715  * \param[in]  rpost  Reference position type to use (NULL = default).
716  * \param[in]  scanner Scanner data structure.
717  * \returns    The created selection element.
718  *
719  * This function handles the creation of a gmx::SelectionTreeElement object for
720  * selection methods that do not take parameters.
721  */
722 SelectionTreeElementPointer _gmx_sel_init_keyword(gmx_ana_selmethod_t*                 method,
723                                                   gmx::SelectionParserValueListPointer args,
724                                                   const char*                          rpost,
725                                                   yyscan_t                             scanner)
726 {
727     return init_keyword_internal(method, gmx::eStringMatchType_Auto, std::move(args), rpost, scanner);
728 }
729
730 /*!
731  * \param[in]  method    Method to use.
732  * \param[in]  matchType String matching type.
733  * \param[in]  args      Pointer to the first argument.
734  * \param[in]  rpost     Reference position type to use (NULL = default).
735  * \param[in]  scanner   Scanner data structure.
736  * \returns    The created selection element.
737  *
738  * This function handles the creation of a gmx::SelectionTreeElement object for
739  * keyword string matching.
740  */
741 SelectionTreeElementPointer _gmx_sel_init_keyword_strmatch(gmx_ana_selmethod_t*          method,
742                                                            gmx::SelectionStringMatchType matchType,
743                                                            gmx::SelectionParserValueListPointer args,
744                                                            const char* rpost,
745                                                            yyscan_t    scanner)
746 {
747     GMX_RELEASE_ASSERT(method->type == STR_VALUE,
748                        "String keyword method called for a non-string-valued method");
749     GMX_RELEASE_ASSERT(args && !args->empty(),
750                        "String keyword matching method called without any values");
751     return init_keyword_internal(method, matchType, std::move(args), rpost, scanner);
752 }
753
754 /*!
755  * \param[in]  method Method to use for initialization.
756  * \param[in]  group  Selection in which the keyword should be evaluated.
757  * \param[in]  rpost  Reference position type to use (NULL = default).
758  * \param[in]  scanner Scanner data structure.
759  * \returns    The created selection element.
760  *
761  * This function handles the creation of a gmx::SelectionTreeElement object for
762  * expressions like "z of ...".
763  */
764 SelectionTreeElementPointer _gmx_sel_init_keyword_of(gmx_ana_selmethod_t*                    method,
765                                                      const gmx::SelectionTreeElementPointer& group,
766                                                      const char*                             rpost,
767                                                      yyscan_t scanner)
768 {
769     // TODO Provide an error if rpost is provided.
770     GMX_UNUSED_VALUE(rpost);
771     return _gmx_sel_init_keyword_evaluator(method, group, scanner);
772 }
773
774 /*!
775  * \param[in]  method Method to use for initialization.
776  * \param[in]  params Pointer to the first parameter.
777  * \param[in]  rpost  Reference position type to use (NULL = default).
778  * \param[in]  scanner Scanner data structure.
779  * \returns    The created selection element.
780  *
781  * This function handles the creation of a gmx::SelectionTreeElement object for
782  * selection methods that take parameters.
783  *
784  * Part of the behavior of the \c same selection keyword is hardcoded into
785  * this function (or rather, into _gmx_selelem_custom_init_same()) to allow the
786  * use of any keyword in \c "same KEYWORD as" without requiring special
787  * handling somewhere else (or sacrificing the simple syntax).
788  */
789 SelectionTreeElementPointer _gmx_sel_init_method(gmx_ana_selmethod_t*                     method,
790                                                  gmx::SelectionParserParameterListPointer params,
791                                                  const char*                              rpost,
792                                                  yyscan_t                                 scanner)
793 {
794     gmx_ana_selcollection_t* sc = _gmx_sel_lexer_selcollection(scanner);
795
796     _gmx_sel_finish_method(scanner);
797     /* The "same" keyword needs some custom massaging of the parameters. */
798     _gmx_selelem_custom_init_same(&method, params, scanner);
799     SelectionTreeElementPointer root(
800             new SelectionTreeElement(SEL_EXPRESSION, _gmx_sel_lexer_get_current_location(scanner)));
801     _gmx_selelem_set_method(root, method, scanner);
802     /* Process the parameters */
803     _gmx_sel_parse_params(*params, root->u.expr.method->nparams, root->u.expr.method->param, root, scanner);
804     set_refpos_type(&sc->pcc, root, rpost);
805
806     return root;
807 }
808
809 /*!
810  * \param[in]  method Modifier to use for initialization.
811  * \param[in]  params Pointer to the first parameter.
812  * \param[in]  sel    Selection element that the modifier should act on.
813  * \param[in]  scanner Scanner data structure.
814  * \returns    The created selection element.
815  *
816  * This function handles the creation of a gmx::SelectionTreeElement object for
817  * selection modifiers.
818  */
819 SelectionTreeElementPointer _gmx_sel_init_modifier(gmx_ana_selmethod_t*                     method,
820                                                    gmx::SelectionParserParameterListPointer params,
821                                                    const gmx::SelectionTreeElementPointer&  sel,
822                                                    yyscan_t                                 scanner)
823 {
824     _gmx_sel_finish_method(scanner);
825     SelectionTreeElementPointer modifier(
826             new SelectionTreeElement(SEL_MODIFIER, _gmx_sel_lexer_get_current_location(scanner)));
827     _gmx_selelem_set_method(modifier, method, scanner);
828     SelectionTreeElementPointer root;
829     if (method->type == NO_VALUE)
830     {
831         SelectionTreeElementPointer child = sel;
832         while (child->next)
833         {
834             child = child->next;
835         }
836         child->next = modifier;
837         root        = sel;
838     }
839     else
840     {
841         params->push_front(SelectionParserParameter::createFromExpression(nullptr, sel));
842         root = modifier;
843     }
844     /* Process the parameters */
845     _gmx_sel_parse_params(*params, modifier->u.expr.method->nparams, modifier->u.expr.method->param,
846                           modifier, scanner);
847
848     return root;
849 }
850
851 /*!
852  * \param[in]  expr    Input selection element for the position calculation.
853  * \param[in]  type    Reference position type or NULL for default.
854  * \param[in]  scanner Scanner data structure.
855  * \returns    The created selection element.
856  *
857  * This function handles the creation of a gmx::SelectionTreeElement object for
858  * evaluation of reference positions.
859  */
860 SelectionTreeElementPointer _gmx_sel_init_position(const gmx::SelectionTreeElementPointer& expr,
861                                                    const char*                             type,
862                                                    yyscan_t                                scanner)
863 {
864     SelectionTreeElementPointer root(
865             new SelectionTreeElement(SEL_EXPRESSION, _gmx_sel_lexer_get_current_location(scanner)));
866     _gmx_selelem_set_method(root, &sm_keyword_pos, scanner);
867     _gmx_selelem_set_kwpos_type(root.get(), type);
868     /* Create the parameters for the parameter parser. */
869     SelectionParserParameterList params;
870     params.push_back(SelectionParserParameter::createFromExpression(nullptr, expr));
871     /* Parse the parameters. */
872     _gmx_sel_parse_params(params, root->u.expr.method->nparams, root->u.expr.method->param, root, scanner);
873
874     return root;
875 }
876
877 /*!
878  * \param[in] x,y,z   Coordinates for the position.
879  * \param[in] scanner Scanner data structure.
880  * \returns   The creates selection element.
881  */
882 SelectionTreeElementPointer _gmx_sel_init_const_position(real x, real y, real z, yyscan_t scanner)
883 {
884     rvec pos{ x, y, z };
885
886     SelectionTreeElementPointer sel(
887             new SelectionTreeElement(SEL_CONST, _gmx_sel_lexer_get_current_location(scanner)));
888     _gmx_selelem_set_vtype(sel, POS_VALUE);
889     _gmx_selvalue_reserve(&sel->v, 1);
890     gmx_ana_pos_init_const(sel->v.u.p, pos);
891     return sel;
892 }
893
894 /*!
895  * \param[in] name  Name of an index group to search for.
896  * \param[in] scanner Scanner data structure.
897  * \returns   The created selection element.
898  *
899  * See gmx_ana_indexgrps_find() for information on how \p name is matched
900  * against the index groups.
901  */
902 SelectionTreeElementPointer _gmx_sel_init_group_by_name(const char* name, yyscan_t scanner)
903 {
904
905     SelectionTreeElementPointer sel(
906             new SelectionTreeElement(SEL_GROUPREF, _gmx_sel_lexer_get_current_location(scanner)));
907     _gmx_selelem_set_vtype(sel, GROUP_VALUE);
908     sel->setName(gmx::formatString("group \"%s\"", name));
909     sel->u.gref.name = gmx_strdup(name);
910     sel->u.gref.id   = -1;
911
912     if (_gmx_sel_lexer_has_groups_set(scanner))
913     {
914         gmx_ana_indexgrps_t*     grps = _gmx_sel_lexer_indexgrps(scanner);
915         gmx_ana_selcollection_t* sc   = _gmx_sel_lexer_selcollection(scanner);
916         sel->resolveIndexGroupReference(grps, sc->gall.isize);
917     }
918
919     return sel;
920 }
921
922 /*!
923  * \param[in] id    Zero-based index number of the group to extract.
924  * \param[in] scanner Scanner data structure.
925  * \returns   The created selection element.
926  */
927 SelectionTreeElementPointer _gmx_sel_init_group_by_id(int id, yyscan_t scanner)
928 {
929     SelectionTreeElementPointer sel(
930             new SelectionTreeElement(SEL_GROUPREF, _gmx_sel_lexer_get_current_location(scanner)));
931     _gmx_selelem_set_vtype(sel, GROUP_VALUE);
932     sel->setName(gmx::formatString("group %d", id));
933     sel->u.gref.name = nullptr;
934     sel->u.gref.id   = id;
935
936     if (_gmx_sel_lexer_has_groups_set(scanner))
937     {
938         gmx_ana_indexgrps_t*     grps = _gmx_sel_lexer_indexgrps(scanner);
939         gmx_ana_selcollection_t* sc   = _gmx_sel_lexer_selcollection(scanner);
940         sel->resolveIndexGroupReference(grps, sc->gall.isize);
941     }
942
943     return sel;
944 }
945
946 /*!
947  * \param[in,out] sel      Value of the variable.
948  * \param         scanner  Scanner data structure.
949  * \returns       The created selection element that references \p sel.
950  *
951  * The reference count of \p sel is updated, but no other modifications are
952  * made.
953  */
954 SelectionTreeElementPointer _gmx_sel_init_variable_ref(const gmx::SelectionTreeElementPointer& sel,
955                                                        yyscan_t scanner)
956 {
957     SelectionTreeElementPointer ref;
958
959     if (sel->v.type == POS_VALUE && sel->type == SEL_CONST)
960     {
961         ref = sel;
962     }
963     else
964     {
965         ref = std::make_shared<SelectionTreeElement>(SEL_SUBEXPRREF,
966                                                      _gmx_sel_lexer_get_current_location(scanner));
967         _gmx_selelem_set_vtype(ref, sel->v.type);
968         ref->setName(sel->name());
969         ref->child = sel;
970     }
971     return ref;
972 }
973
974 /*!
975  * \param[in]  name     Name for the selection
976  *     (if NULL, a default name is constructed).
977  * \param[in]  sel      The selection element that evaluates the selection.
978  * \param      scanner  Scanner data structure.
979  * \returns    The created root selection element.
980  *
981  * This function handles the creation of root (\ref SEL_ROOT)
982  * gmx::SelectionTreeElement objects for selections.
983  */
984 SelectionTreeElementPointer _gmx_sel_init_selection(const char*                             name,
985                                                     const gmx::SelectionTreeElementPointer& sel,
986                                                     yyscan_t                                scanner)
987 {
988     if (sel->v.type != POS_VALUE)
989     {
990         /* FIXME: Better handling of this error */
991         GMX_THROW(gmx::InternalError("Each selection must evaluate to a position"));
992     }
993
994     SelectionTreeElementPointer root(
995             new SelectionTreeElement(SEL_ROOT, _gmx_sel_lexer_get_current_location(scanner)));
996     root->child = sel;
997     if (name)
998     {
999         root->setName(name);
1000     }
1001     /* Update the flags */
1002     _gmx_selelem_update_flags(root);
1003     gmx::ExceptionInitializer errors("Invalid index group reference(s)");
1004     root->checkUnsortedAtoms(true, &errors);
1005     if (errors.hasNestedExceptions())
1006     {
1007         GMX_THROW(gmx::InconsistentInputError(errors));
1008     }
1009
1010     root->fillNameIfMissing(_gmx_sel_lexer_pselstr(scanner));
1011
1012     /* Print out some information if the parser is interactive */
1013     gmx::TextWriter* statusWriter = _gmx_sel_lexer_get_status_writer(scanner);
1014     if (statusWriter != nullptr)
1015     {
1016         const std::string message =
1017                 gmx::formatString("Selection '%s' parsed", _gmx_sel_lexer_pselstr(scanner));
1018         statusWriter->writeLine(message);
1019     }
1020
1021     return root;
1022 }
1023
1024
1025 /*!
1026  * \param[in]  name     Name of the variable.
1027  * \param[in]  expr     The selection element that evaluates the variable.
1028  * \param      scanner  Scanner data structure.
1029  * \returns    The created root selection element.
1030  *
1031  * This function handles the creation of root gmx::SelectionTreeElement objects
1032  * for variable assignments. A \ref SEL_ROOT element and a \ref SEL_SUBEXPR
1033  * element are both created.
1034  */
1035 SelectionTreeElementPointer _gmx_sel_assign_variable(const char*                             name,
1036                                                      const gmx::SelectionTreeElementPointer& expr,
1037                                                      yyscan_t scanner)
1038 {
1039     gmx_ana_selcollection_t*    sc      = _gmx_sel_lexer_selcollection(scanner);
1040     const char*                 pselstr = _gmx_sel_lexer_pselstr(scanner);
1041     SelectionTreeElementPointer root;
1042
1043     _gmx_selelem_update_flags(expr);
1044     /* Check if this is a constant non-group value */
1045     if (expr->type == SEL_CONST && expr->v.type != GROUP_VALUE)
1046     {
1047         /* If so, just assign the constant value to the variable */
1048         sc->symtab->addVariable(name, expr);
1049     }
1050     /* Check if we are assigning a variable to another variable */
1051     else if (expr->type == SEL_SUBEXPRREF)
1052     {
1053         /* If so, make a simple alias */
1054         sc->symtab->addVariable(name, expr->child);
1055     }
1056     else
1057     {
1058         SelectionLocation location(_gmx_sel_lexer_get_current_location(scanner));
1059         /* Create the root element */
1060         root = std::make_shared<SelectionTreeElement>(SEL_ROOT, location);
1061         root->setName(name);
1062         /* Create the subexpression element */
1063         root->child = std::make_shared<SelectionTreeElement>(SEL_SUBEXPR, location);
1064         root->child->setName(name);
1065         _gmx_selelem_set_vtype(root->child, expr->v.type);
1066         root->child->child = expr;
1067         /* Update flags */
1068         _gmx_selelem_update_flags(root);
1069         gmx::ExceptionInitializer errors("Invalid index group reference(s)");
1070         root->checkUnsortedAtoms(true, &errors);
1071         if (errors.hasNestedExceptions())
1072         {
1073             GMX_THROW(gmx::InconsistentInputError(errors));
1074         }
1075         /* Add the variable to the symbol table */
1076         sc->symtab->addVariable(name, root->child);
1077     }
1078     srenew(sc->varstrs, sc->nvars + 1);
1079     sc->varstrs[sc->nvars] = gmx_strdup(pselstr);
1080     ++sc->nvars;
1081     gmx::TextWriter* statusWriter = _gmx_sel_lexer_get_status_writer(scanner);
1082     if (statusWriter != nullptr)
1083     {
1084         const std::string message = gmx::formatString("Variable '%s' parsed", pselstr);
1085         statusWriter->writeLine(message);
1086     }
1087     return root;
1088 }
1089
1090 /*!
1091  * \param         sel   Selection to append (can be NULL, in which
1092  *   case nothing is done).
1093  * \param         last  Last selection, or NULL if not present or not known.
1094  * \param         scanner  Scanner data structure.
1095  * \returns       The last selection after the append.
1096  *
1097  * Appends \p sel after the last root element, and returns either \p sel
1098  * (if it was non-NULL) or the last element (if \p sel was NULL).
1099  */
1100 SelectionTreeElementPointer _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(sel.get(), _gmx_sel_lexer_pselstr(scanner)));
1136             sc->sel.push_back(std::move(selPtr));
1137         }
1138     }
1139     /* Clear the selection string now that we've saved it */
1140     _gmx_sel_lexer_clear_pselstr(scanner);
1141     return last;
1142 }
1143
1144 /*!
1145  * \param[in] scanner Scanner data structure.
1146  * \returns   true if the parser should finish, false if parsing should
1147  *   continue.
1148  *
1149  * This function is called always after _gmx_sel_append_selection() to
1150  * check whether a sufficient number of selections has already been provided.
1151  * This is used to terminate interactive parsers when the correct number of
1152  * selections has been provided.
1153  */
1154 bool _gmx_sel_parser_should_finish(yyscan_t scanner)
1155 {
1156     gmx_ana_selcollection_t* sc = _gmx_sel_lexer_selcollection(scanner);
1157     return gmx::ssize(sc->sel) == _gmx_sel_lexer_exp_selcount(scanner);
1158 }