Replace compat::make_unique with std::make_unique
[alexxy/gromacs.git] / src / gromacs / selection / parsetree.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2009,2010,2011,2012,2013,2014,2015,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  * Handling of intermediate selection parser data.
38  *
39  * The data types declared in this header are used by the parser to store
40  * intermediate data when constructing method expressions.
41  * In particular, the parameters for the method are stored.
42  * The intermediate data is freed once a gmx::SelectionTreeElement object can
43  * be constructed.
44  *
45  * This is an implementation header: there should be no need to use it outside
46  * this directory.
47  *
48  * \author Teemu Murtola <teemu.murtola@gmail.com>
49  * \ingroup module_selection
50  */
51 #ifndef GMX_SELECTION_PARSETREE_H
52 #define GMX_SELECTION_PARSETREE_H
53
54 #include <exception>
55 #include <list>
56 #include <memory>
57 #include <string>
58
59 #include "gromacs/math/vec.h"
60 #include "gromacs/math/vectypes.h"
61 #include "gromacs/utility/gmxassert.h"
62 #include "gromacs/utility/real.h"
63
64 #include "selelem.h"
65 #include "selvalue.h"
66
67 struct gmx_ana_indexgrps_t;
68 struct gmx_ana_selmethod_t;
69 struct gmx_ana_selparam_t;
70
71 namespace gmx
72 {
73
74 //! \cond internal
75 /*! \internal \brief
76  * String matching mode for string keyword expressions.
77  *
78  * \ingroup module_selection
79  */
80 enum SelectionStringMatchType
81 {
82     eStringMatchType_Auto,              //!< Deduce from the string.
83     eStringMatchType_Exact,             //!< Match as a literal string.
84     eStringMatchType_Wildcard,          //!< Match using ? and * as wildcards.
85     eStringMatchType_RegularExpression  //!< Match using regular expressions.
86 };
87 /*! \endcond */
88
89 class SelectionParserValue;
90
91 //! Container for a list of SelectionParserValue objects.
92 typedef std::list<SelectionParserValue>
93     SelectionParserValueList;
94 //! Smart pointer type for managing a SelectionParserValueList.
95 typedef std::unique_ptr<SelectionParserValueList>
96     SelectionParserValueListPointer;
97
98 /*! \internal
99  * \brief
100  * Describes a parsed value, possibly resulting from expression evaluation.
101  *
102  * All factory methods and the constructors may throw an std::bad_alloc if
103  * out of memory.
104  *
105  * \ingroup module_selection
106  */
107 class SelectionParserValue
108 {
109     public:
110         //! Allocates and initializes an empty value list.
111         static SelectionParserValueListPointer createList()
112         {
113             return std::make_unique<SelectionParserValueList>();
114         }
115         /*! \brief
116          * Allocates and initializes a value list with a single value.
117          *
118          * \param[in] value  Initial value to put in the list.
119          * \returns   Pointer to a new value list that contains \p value.
120          */
121         static SelectionParserValueListPointer
122         createList(const SelectionParserValue &value)
123         {
124             SelectionParserValueListPointer list(new SelectionParserValueList);
125             list->push_back(value);
126             return list;
127         }
128         /*! \brief
129          * Allocates and initializes an expression value.
130          *
131          * \param[in] expr  Root of the expression tree to assign to the value.
132          * \returns   The newly created value.
133          */
134         static SelectionParserValue
135         createExpr(const gmx::SelectionTreeElementPointer &expr)
136         {
137             return SelectionParserValue(expr);
138         }
139         /*! \brief
140          * Allocates and initializes a constant integer value.
141          *
142          * \param[in] value    Integer value to assign to the value.
143          * \param[in] location Location of the value.
144          * \returns   The newly created value.
145          */
146         static SelectionParserValue
147         createInteger(int value, const SelectionLocation &location)
148         {
149             SelectionParserValue result(INT_VALUE, location);
150             result.u.i.i1 = result.u.i.i2 = value;
151             return result;
152         }
153         /*! \brief
154          * Allocates and initializes a constant integer range value.
155          *
156          * \param[in] from     Beginning of the range to assign to the value.
157          * \param[in] to       End of the range to assign to the value.
158          * \param[in] location Location of the value.
159          * \returns   The newly created value.
160          */
161         static SelectionParserValue
162         createIntegerRange(int from, int to, const SelectionLocation &location)
163         {
164             SelectionParserValue result(INT_VALUE, location);
165             result.u.i.i1 = from;
166             result.u.i.i2 = to;
167             return result;
168         }
169         /*! \brief
170          * Allocates and initializes a constant floating-point value.
171          *
172          * \param[in] value    Floating-point value to assign to the value.
173          * \param[in] location Location of the value.
174          * \returns   The newly created value.
175          */
176         static SelectionParserValue
177         createReal(real value, const SelectionLocation &location)
178         {
179             SelectionParserValue result(REAL_VALUE, location);
180             result.u.r.r1 = result.u.r.r2 = value;
181             return result;
182         }
183         /*! \brief
184          * Allocates and initializes a constant floating-point range value.
185          *
186          * \param[in] from     Beginning of the range to assign to the value.
187          * \param[in] to       End of the range to assign to the value.
188          * \param[in] location Location of the value.
189          * \returns   The newly created value.
190          */
191         static SelectionParserValue
192         createRealRange(real from, real to, const SelectionLocation &location)
193         {
194             SelectionParserValue result(REAL_VALUE, location);
195             result.u.r.r1 = from;
196             result.u.r.r2 = to;
197             return result;
198         }
199         /*! \brief
200          * Allocates and initializes a constant string value.
201          *
202          * \param[in] value    String to assign to the value.
203          * \param[in] location Location of the value.
204          * \returns   The newly created value.
205          */
206         static SelectionParserValue
207         createString(const char *value, const SelectionLocation &location)
208         {
209             SelectionParserValue result(STR_VALUE, location);
210             result.str = value;
211             return result;
212         }
213         /*! \brief
214          * Allocates and initializes a constant position value.
215          *
216          * \param[in] value    Position vector to assign to the value.
217          * \param[in] location Location of the value.
218          * \returns   The newly created value.
219          */
220         static SelectionParserValue
221         createPosition(rvec value, const SelectionLocation &location)
222         {
223             SelectionParserValue result(POS_VALUE, location);
224             copy_rvec(value, result.u.x);
225             return result;
226         }
227
228         //! Returns the location of this value in the parsed selection text.
229         const SelectionLocation &location() const { return location_; }
230         //! Returns true if the value comes from expression evaluation.
231         bool hasExpressionValue() const { return static_cast<bool>(expr); }
232
233         //! Returns the string value (\a type must be ::STR_VALUE).
234         const std::string &stringValue() const
235         {
236             GMX_ASSERT(type == STR_VALUE && !hasExpressionValue(),
237                        "Attempted to retrieve string value from a non-string value");
238             return str;
239         }
240
241         // TODO: boost::any or similar could be nicer for the implementation.
242         //! Type of the value.
243         e_selvalue_t                     type;
244         //! Expression pointer if the value is the result of an expression.
245         gmx::SelectionTreeElementPointer expr;
246         //! String value for \a type ::STR_VALUE.
247         std::string                      str;
248         //! The actual value if \a expr is NULL and \a type is not ::STR_VALUE.
249         union {
250             //! The integer value/range (\a type ::INT_VALUE).
251             struct {
252                 //! Beginning of the range.
253                 int             i1;
254                 //! End of the range; equals \a i1 for a single integer.
255                 int             i2;
256             }                   i;
257             //! The real value/range (\a type ::REAL_VALUE).
258             struct {
259                 //! Beginning of the range.
260                 real            r1;
261                 //! End of the range; equals \a r1 for a single number.
262                 real            r2;
263             }                   r;
264             //! The position value (\a type ::POS_VALUE).
265             rvec                x;
266         }                       u;
267
268     private:
269         /*! \brief
270          * Initializes a new value.
271          *
272          * \param[in] type     Type for the new value.
273          * \param[in] location Location for the value.
274          */
275         SelectionParserValue(e_selvalue_t type, const SelectionLocation &location);
276         /*! \brief
277          * Initializes a new expression value.
278          *
279          * \param[in] expr  Expression for the value.
280          */
281         explicit SelectionParserValue(const gmx::SelectionTreeElementPointer &expr);
282
283         //! Location of the value in the parsed text.
284         SelectionLocation       location_;
285 };
286
287 class SelectionParserParameter;
288
289 //! Container for a list of SelectionParserParameter objects.
290 typedef std::list<SelectionParserParameter>
291     SelectionParserParameterList;
292 //! Smart pointer type for managing a SelectionParserParameterList.
293 typedef std::unique_ptr<SelectionParserParameterList>
294     SelectionParserParameterListPointer;
295
296 /*! \internal \brief
297  * Describes a parsed method parameter.
298  *
299  * \ingroup module_selection
300  */
301 class SelectionParserParameter
302 {
303     public:
304         // Default move constructor and assignment. Only needed for old compilers.
305         //! \cond
306         SelectionParserParameter(SelectionParserParameter &&o) noexcept
307             : name_(std::move(o.name_)), location_(o.location_),
308               values_(std::move(o.values_))
309         {
310         }
311
312         SelectionParserParameter &operator=(SelectionParserParameter &&o) noexcept
313         {
314             name_     = std::move(o.name_);
315             location_ = o.location_;
316             values_   = std::move(o.values_);
317             return *this;
318         }
319         //! \endcond
320
321         //! Allocates and initializes an empty parameter list.
322         static SelectionParserParameterListPointer createList()
323         {
324             return std::make_unique<SelectionParserParameterList>();
325         }
326         /*! \brief
327          * Allocates and initializes a parsed method parameter.
328          *
329          * \param[in] name     Name for the new parameter (can be NULL).
330          * \param[in] values   List of values for the parameter.
331          * \param[in] location Location of the parameter.
332          * \returns   Pointer to the newly allocated parameter.
333          * \throws    std::bad_alloc if out of memory.
334          */
335         static SelectionParserParameter
336         create(const char *name, SelectionParserValueListPointer values,
337                const SelectionLocation &location)
338         {
339             return SelectionParserParameter(name, std::move(values), location);
340         }
341         //! \copydoc create(const char *, SelectionParserValueListPointer, const SelectionLocation &)
342         static SelectionParserParameter
343         create(const std::string &name, SelectionParserValueListPointer values,
344                const SelectionLocation &location)
345         {
346             return SelectionParserParameter(name.c_str(), std::move(values), location);
347         }
348         /*! \brief
349          * Allocates and initializes a parsed method parameter.
350          *
351          * \param[in] name     Name for the new parameter (can be NULL).
352          * \param[in] value    Value for the parameter.
353          * \param[in] location Location of the parameter.
354          * \returns   Pointer to the newly allocated parameter.
355          * \throws    std::bad_alloc if out of memory.
356          *
357          * This overload is a convenience wrapper for the case when creating
358          * parameters outside the actual Bison parser and only a single value
359          * is necessary.
360          */
361         static SelectionParserParameter
362         create(const char *name, const SelectionParserValue &value,
363                const SelectionLocation &location)
364         {
365             return create(name, SelectionParserValue::createList(value), location);
366         }
367         /*! \brief
368          * Allocates and initializes a parsed method parameter.
369          *
370          * \param[in] name    Name for the new parameter (can be NULL).
371          * \param[in] expr    Expression value for the parameter.
372          * \returns   Pointer to the newly allocated parameter.
373          * \throws    std::bad_alloc if out of memory.
374          *
375          * This overload is a convenience wrapper for the case when creating
376          * parameters outside the actual Bison parser and only a single
377          * expression value is necessary.
378          */
379         static SelectionParserParameter
380         createFromExpression(const char                        *name,
381                              const SelectionTreeElementPointer &expr)
382         {
383             return create(name, SelectionParserValue::createExpr(expr),
384                           expr->location());
385         }
386         //! \copydoc createFromExpression(const char *, const SelectionTreeElementPointer &)
387         static SelectionParserParameter
388         createFromExpression(const std::string                 &name,
389                              const SelectionTreeElementPointer &expr)
390         {
391             return create(name.c_str(), SelectionParserValue::createExpr(expr),
392                           expr->location());
393         }
394
395         //! Returns the name of the parameter (may be empty).
396         const std::string &name() const { return name_; }
397         //! Returns the location of this parameter in the parsed selection text.
398         const SelectionLocation        &location() const { return location_; }
399         //! Returns the values for the parameter.
400         const SelectionParserValueList &values() const { return *values_; }
401
402     private:
403         /*! \brief
404          * Initializes a parsed method parameter.
405          *
406          * \param[in] name     Name for the new parameter (can be NULL).
407          * \param[in] values   List of values for the parameter.
408          * \param[in] location Location of the parameter.
409          * \throws    std::bad_alloc if out of memory.
410          */
411         SelectionParserParameter(const char                      *name,
412                                  SelectionParserValueListPointer  values,
413                                  const SelectionLocation         &location);
414
415         //! Name of the parameter.
416         std::string                     name_;
417         //! Location of the parameter in the parsed text.
418         SelectionLocation               location_;
419
420         // TODO: Make private, there is only one direct user.
421     public:
422         //! Values for this parameter.
423         SelectionParserValueListPointer values_;
424 };
425
426 } // namespace gmx
427
428 /*! \brief
429  * Handles exceptions caught within the Bison code.
430  *
431  * \retval `true`  if the parser should attempt error recovery.
432  * \retval `false` if the parser should immediately abort.
433  *
434  * This function is called whenever an exception is caught within Bison
435  * actions.  Since exceptions cannot propagate through Bison code, the
436  * exception is saved (potentially with some extra context information) so that
437  * the caller of the parser can rethrow the exception.
438  *
439  * If it is possible to recover from the exception, then the function returns
440  * `true`, and Bison enters error recovery state.  At the end of the recovery,
441  * _gmx_selparser_handle_error() is called.
442  * If this function returns false, then Bison immediately aborts the parsing
443  * so that the caller can rethrow the exception.
444  */
445 bool
446 _gmx_selparser_handle_exception(void *scanner, std::exception *ex);
447 /*! \brief
448  * Handles errors in the selection parser.
449  *
450  * \returns `true` if parsing can continue with the next selection.
451  * \throws  std::bad_alloc if out of memory during the error processing.
452  * \throws  unspecified    Can throw the stored exception if recovery from that
453  *     exception is not possible.
454  *
455  * This function is called during error recovery, after Bison has discarded all
456  * the symbols for the erroneous selection.
457  * At this point, the full selection that caused the error is known, and can be
458  * added to the error context.
459  *
460  * For an interactive parser, this function returns `true` to let the parsing
461  * continue with the next selection, or to let the user enter the next
462  * selection, if it was possible to recover from the exception.
463  * For other cases, this will either rethrow the original exception with added
464  * context, or return `false` after adding the context to the error reporter.
465  * Any exceptions thrown from this method are again caught by Bison and result
466  * in termination of the parsing; the caller can then rethrow them.
467  */
468 bool
469 _gmx_selparser_handle_error(void *scanner);
470
471 /** Propagates the flags for selection elements. */
472 void
473 _gmx_selelem_update_flags(const gmx::SelectionTreeElementPointer &sel);
474
475 /** Initializes the method parameter data of \ref SEL_EXPRESSION and
476  * \ref SEL_MODIFIER elements. */
477 void
478 _gmx_selelem_init_method_params(const gmx::SelectionTreeElementPointer &sel,
479                                 void                                   *scanner);
480 /** Initializes the method for a \ref SEL_EXPRESSION selection element. */
481 void
482 _gmx_selelem_set_method(const gmx::SelectionTreeElementPointer &sel,
483                         struct gmx_ana_selmethod_t *method, void *scanner);
484
485 /* An opaque pointer. */
486 #ifndef YY_TYPEDEF_YY_SCANNER_T
487 #define YY_TYPEDEF_YY_SCANNER_T
488 typedef void* yyscan_t;
489 #endif
490 /** \brief Creates a gmx::SelectionTreeElement for arithmetic expression evaluation.
491  *
492  * \param[in]  left    Selection element for the left hand side.
493  * \param[in]  right   Selection element for the right hand side.
494  * \param[in]  op      String representation of the operator.
495  * \param[in]  scanner Scanner data structure.
496  * \returns    The created selection element.
497  *
498  * This function handles the creation of a gmx::SelectionTreeElement object for
499  * arithmetic expressions.
500  */
501 gmx::SelectionTreeElementPointer
502 _gmx_sel_init_arithmetic(const gmx::SelectionTreeElementPointer &left,
503                          const gmx::SelectionTreeElementPointer &right,
504                          char op, yyscan_t scanner);
505 /** Creates a gmx::SelectionTreeElement for comparsion expression evaluation. */
506 gmx::SelectionTreeElementPointer
507 _gmx_sel_init_comparison(const gmx::SelectionTreeElementPointer &left,
508                          const gmx::SelectionTreeElementPointer &right,
509                          const char *cmpop, void *scanner);
510 /** Creates a gmx::SelectionTreeElement for a keyword expression from the parsed data. */
511 gmx::SelectionTreeElementPointer
512 _gmx_sel_init_keyword(struct gmx_ana_selmethod_t *method,
513                       gmx::SelectionParserValueListPointer args,
514                       const char *rpost, void *scanner);
515 /** Creates a gmx::SelectionTreeElement for string-matching keyword expression. */
516 gmx::SelectionTreeElementPointer
517 _gmx_sel_init_keyword_strmatch(struct gmx_ana_selmethod_t *method,
518                                gmx::SelectionStringMatchType matchType,
519                                gmx::SelectionParserValueListPointer args,
520                                const char *rpost, void *scanner);
521 /** Creates a gmx::SelectionTreeElement for "keyword of" expression. */
522 gmx::SelectionTreeElementPointer
523 _gmx_sel_init_keyword_of(struct gmx_ana_selmethod_t *method,
524                          const gmx::SelectionTreeElementPointer &group,
525                          const char *rpost, void *scanner);
526 /** Creates a gmx::SelectionTreeElement for a method expression from the parsed data. */
527 gmx::SelectionTreeElementPointer
528 _gmx_sel_init_method(struct gmx_ana_selmethod_t *method,
529                      gmx::SelectionParserParameterListPointer params,
530                      const char *rpost, void *scanner);
531 /** Creates a gmx::SelectionTreeElement for a modifier expression from the parsed data. */
532 gmx::SelectionTreeElementPointer
533 _gmx_sel_init_modifier(struct gmx_ana_selmethod_t              *mod,
534                        gmx::SelectionParserParameterListPointer params,
535                        const gmx::SelectionTreeElementPointer  &sel,
536                        void                                    *scanner);
537 /** Creates a gmx::SelectionTreeElement for evaluation of reference positions. */
538 gmx::SelectionTreeElementPointer
539 _gmx_sel_init_position(const gmx::SelectionTreeElementPointer &expr,
540                        const char *type, void *scanner);
541
542 /** Creates a gmx::SelectionTreeElement for a constant position. */
543 gmx::SelectionTreeElementPointer
544 _gmx_sel_init_const_position(real x, real y, real z, void *scanner);
545 /** Creates a gmx::SelectionTreeElement for a index group expression using group name. */
546 gmx::SelectionTreeElementPointer
547 _gmx_sel_init_group_by_name(const char *name, void *scanner);
548 /** Creates a gmx::SelectionTreeElement for a index group expression using group index. */
549 gmx::SelectionTreeElementPointer
550 _gmx_sel_init_group_by_id(int id, void *scanner);
551 /** Creates a gmx::SelectionTreeElement for a variable reference */
552 gmx::SelectionTreeElementPointer
553 _gmx_sel_init_variable_ref(const gmx::SelectionTreeElementPointer &sel,
554                            void                                   *scanner);
555
556 /** Creates a root gmx::SelectionTreeElement for a selection. */
557 gmx::SelectionTreeElementPointer
558 _gmx_sel_init_selection(const char                             *name,
559                         const gmx::SelectionTreeElementPointer &sel,
560                         void                                   *scanner);
561 /** Creates a root gmx::SelectionTreeElement elements for a variable assignment. */
562 gmx::SelectionTreeElementPointer
563 _gmx_sel_assign_variable(const char                             *name,
564                          const gmx::SelectionTreeElementPointer &expr,
565                          void                                   *scanner);
566 /** Appends a root gmx::SelectionTreeElement to a selection collection. */
567 gmx::SelectionTreeElementPointer
568 _gmx_sel_append_selection(const gmx::SelectionTreeElementPointer &sel,
569                           gmx::SelectionTreeElementPointer        last,
570                           void                                   *scanner);
571 /** Check whether the parser should finish. */
572 bool
573 _gmx_sel_parser_should_finish(void *scanner);
574
575 /* In params.c */
576 /** Initializes an array of parameters based on input from the selection parser. */
577 void
578 _gmx_sel_parse_params(const gmx::SelectionParserParameterList &params,
579                       int nparam, struct gmx_ana_selparam_t *param,
580                       const gmx::SelectionTreeElementPointer &root,
581                       void *scanner);
582
583 #endif