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