2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2009-2017, The GROMACS development team.
5 * Copyright (c) 2019, by the GROMACS development team, led by
6 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7 * and including many others, as listed in the AUTHORS file in the
8 * top-level source directory and at http://www.gromacs.org.
10 * GROMACS is free software; you can redistribute it and/or
11 * modify it under the terms of the GNU Lesser General Public License
12 * as published by the Free Software Foundation; either version 2.1
13 * of the License, or (at your option) any later version.
15 * GROMACS is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18 * Lesser General Public License for more details.
20 * You should have received a copy of the GNU Lesser General Public
21 * License along with GROMACS; if not, see
22 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
25 * If you want to redistribute modifications to GROMACS, please
26 * consider that scientific software is very special. Version
27 * control is crucial - bugs must be traceable. We will be happy to
28 * consider code for inclusion in the official distribution, but
29 * derived work must not be called official GROMACS. Details are found
30 * in the README & COPYING files - if they are missing, get the
31 * official version at http://www.gromacs.org.
33 * To help us fund GROMACS development, we humbly ask that you cite
34 * the research papers on the package. Check out http://www.gromacs.org.
38 * Declares gmx::SelectionTreeElement and related things.
40 * The selection element trees constructed by the parser and the compiler
41 * are described on the respective pages:
42 * \ref page_module_selection_parser and \ref page_module_selection_compiler.
44 * This is an implementation header: there should be no need to use it outside
47 * \author Teemu Murtola <teemu.murtola@gmail.com>
48 * \ingroup module_selection
50 #ifndef GMX_SELECTION_SELELEM_H
51 #define GMX_SELECTION_SELELEM_H
56 #include "gromacs/selection/indexutil.h"
57 #include "gromacs/utility/classhelpers.h"
58 #include "gromacs/utility/real.h"
62 struct gmx_ana_poscalc_t
;
63 struct gmx_ana_selparam_t
;
64 struct gmx_ana_selmethod_t
;
66 struct gmx_sel_evaluate_t
;
67 struct gmx_sel_mempool_t
;
69 struct t_compiler_data
;
73 class SelectionTreeElement
;
74 struct SelectionTopologyProperties
;
76 //! Smart pointer type for selection tree element pointers.
77 typedef std::shared_ptr
<SelectionTreeElement
> SelectionTreeElementPointer
;
80 /********************************************************************/
81 /*! \name Enumerations for expression types
82 ********************************************************************/
85 /** Defines the type of a gmx::SelectionTreeElement object. */
88 /** Constant-valued expression. */
90 /** Method expression that requires evaluation. */
92 /** Boolean expression. */
94 /** Arithmetic expression. */
96 /** Root node of the evaluation tree. */
98 /** Subexpression that may be referenced several times. */
100 /** Reference to a subexpression. */
102 /** Unresolved reference to an external group. */
104 /** Post-processing of selection value. */
108 /** Defines the boolean operation of gmx::SelectionTreeElement objects with type \ref SEL_BOOLEAN. */
111 BOOL_NOT
, /**< Not */
112 BOOL_AND
, /**< And */
114 BOOL_XOR
/**< Xor (not implemented). */
117 /** Defines the arithmetic operation of gmx::SelectionTreeElement objects with type \ref SEL_ARITHMETIC. */
120 ARITH_PLUS
, /**< Addition (`+`) */
121 ARITH_MINUS
, /**< Subtraction (`-`) */
122 ARITH_NEG
, /**< Unary `-` */
123 ARITH_MULT
, /**< Multiplication (`*`) */
124 ARITH_DIV
, /**< Division (`/`) */
125 ARITH_EXP
/**< Power (`^`) */
128 /** Returns a string representation of the type of a gmx::SelectionTreeElement. */
129 extern const char* _gmx_selelem_type_str(const gmx::SelectionTreeElement
& sel
);
130 /** Returns a string representation of the boolean type of a \ref SEL_BOOLEAN gmx::SelectionTreeElement. */
131 extern const char* _gmx_selelem_boolean_type_str(const gmx::SelectionTreeElement
& sel
);
132 /** Returns a string representation of the type of a \c gmx_ana_selvalue_t. */
133 extern const char* _gmx_sel_value_type_str(const gmx_ana_selvalue_t
* val
);
138 /********************************************************************/
139 /*! \name Selection expression flags
140 * \anchor selelem_flags
141 ********************************************************************/
144 * Selection value flags are set.
146 * If this flag is set, the flags covered by \ref SEL_VALFLAGMASK
147 * have been set properly for the element.
149 #define SEL_FLAGSSET 1
151 * The element evaluates to a single value.
153 * This flag is always set for \ref GROUP_VALUE elements.
155 #define SEL_SINGLEVAL 2
157 * The element evaluates to one value for each input atom.
159 #define SEL_ATOMVAL 4
161 * The element evaluates to an arbitrary number of values.
163 #define SEL_VARNUMVAL 8
165 * The element (or one of its children) is dynamic.
167 #define SEL_DYNAMIC 16
169 * The element may contain atom indices in an unsorted order.
171 #define SEL_UNSORTED 32
173 * Mask that covers the flags that describe the number of values.
175 #define SEL_VALTYPEMASK (SEL_SINGLEVAL | SEL_ATOMVAL | SEL_VARNUMVAL)
177 * Mask that covers the flags that describe the value type.
179 #define SEL_VALFLAGMASK (SEL_FLAGSSET | SEL_VALTYPEMASK | SEL_DYNAMIC)
181 * Data has been allocated for the \p v.u union.
183 * If not set, the \p v.u.ptr points to data allocated externally.
184 * This is the case if the value of the element is used as a parameter
185 * for a selection method or if the element evaluates the final value of
188 * Even if the flag is set, \p v.u.ptr can be NULL during initialization.
191 * This flag overlaps with the function of \p v.nalloc field, and could
192 * probably be removed, making memory management simpler. Currently, the
193 * \p v.nalloc field is not kept up-to-date in all cases when this flag
194 * is changed and is used in places where this flag is not, so this would
195 * require a careful investigation of the selection code.
197 #define SEL_ALLOCVAL (1 << 8)
199 * Data has been allocated for the group/position structure.
201 * If not set, the memory allocated for fields in \p v.u.g or \p v.u.p is
202 * managed externally.
204 * This field has no effect if the value type is not \ref GROUP_VALUE or
205 * \ref POS_VALUE, but should not be set.
207 #define SEL_ALLOCDATA (1 << 9)
209 * \p method->init_frame should be called for the frame.
211 #define SEL_INITFRAME (1 << 10)
213 * Parameter has been evaluated for the current frame.
215 * This flag is set for children of \ref SEL_EXPRESSION elements (which
216 * describe method parameters) after the element has been evaluated for the
218 * It is not set for \ref SEL_ATOMVAL elements, because they may need to
219 * be evaluated multiple times.
221 #define SEL_EVALFRAME (1 << 11)
223 * \p method->init has been called.
225 #define SEL_METHODINIT (1 << 12)
227 * \p method->outinit has been called.
229 * This flag is also used for \ref SEL_SUBEXPRREF elements.
231 #define SEL_OUTINIT (1 << 13)
238 class ExceptionInitializer
;
242 * Function pointer for evaluating a gmx::SelectionTreeElement.
244 typedef void (*sel_evalfunc
)(struct gmx_sel_evaluate_t
* data
,
245 const SelectionTreeElementPointer
& sel
,
251 * Stores the location of a selection element in the selection text.
253 * The location is stored as a range in the pretty-printed selection text
254 * (where whitespace has been sanitized), and can be used to extract that text
255 * for error messages and other diagnostic purposes.
256 * During parsing, the extraction is done with _gmx_sel_lexer_get_text().
258 * This needs to be a plain C struct for Bison to properly deal with it.
260 struct SelectionLocation
262 //! Returns an empty location.
263 static SelectionLocation
createEmpty()
265 SelectionLocation empty
= { 0, 0 };
269 //! Start index of the string where this element has been parsed from.
271 //! End index of the string where this element has been parsed from.
276 * Represents an element of a selection expression.
278 class SelectionTreeElement
282 * Allocates memory and performs common initialization.
284 * \param[in] type Type of selection element to create.
285 * \param[in] location Location of the element.
287 * \a type is set to \p type,
288 * \a v::type is set to \ref GROUP_VALUE for boolean and comparison
289 * expressions and \ref NO_VALUE for others, and
290 * \ref SEL_ALLOCVAL is set for non-root elements (\ref SEL_ALLOCDATA
291 * is also set for \ref SEL_BOOLEAN elements).
292 * All the pointers are set to NULL.
294 SelectionTreeElement(e_selelem_t type
, const SelectionLocation
& location
);
295 ~SelectionTreeElement();
297 //! Frees the memory allocated for the \a v union.
299 //! Frees the memory allocated for the \a u union.
300 void freeExpressionData();
301 /* In compiler.cpp */
303 * Frees the memory allocated for the selection compiler.
305 * This function only frees the data for the given selection, not its
306 * children. It is safe to call the function when compiler data has
307 * not been allocated or has already been freed; in such a case,
310 void freeCompilerData();
313 * Reserves memory for value from a memory pool.
315 * \param[in] count Number of values to reserve memory for.
317 * Reserves memory for the values of this element from the \a mempool
319 * If no memory pool is set, nothing is done.
321 void mempoolReserve(int count
);
323 * Releases memory pool used for value.
325 * Releases the memory allocated for the values of this element from the
326 * \a mempool memory pool.
327 * If no memory pool is set, nothing is done.
329 void mempoolRelease();
331 //! Returns the name of the element.
332 const std::string
& name() const { return name_
; }
333 //! Returns the location of the element.
334 const SelectionLocation
& location() const { return location_
; }
337 * Sets the name of the element.
339 * \param[in] name Name to set (can be NULL).
340 * \throws std::bad_alloc if out of memory.
342 void setName(const char* name
) { name_
= (name
!= nullptr ? name
: ""); }
343 //! \copydoc setName(const char *)
344 void setName(const std::string
& name
) { name_
= name
; }
346 * Sets the name of a root element if it is missing.
348 * \param[in] selectionText Full selection text to use as a fallback.
349 * \throws std::bad_alloc if out of memory.
351 * If index groups have not yet been set and the selection is a result
352 * of a group reference, the name may still be empty after this call.
354 * Strong exception safety guarantee.
356 void fillNameIfMissing(const char* selectionText
);
359 * Returns which topology properties the selection element subtree requires
362 * \returns List of topology properties required for evaluation.
364 SelectionTopologyProperties
requiredTopologyProperties() const;
366 * Checks that this element and its children do not contain unsupported
367 * elements with unsorted atoms.
369 * \param[in] bUnsortedAllowed Whether this element's parents allow it
370 * to have unsorted atoms.
371 * \param errors Object for reporting any error messages.
372 * \throws std::bad_alloc if out of memory.
374 * Errors are reported as nested exceptions in \p errors.
376 void checkUnsortedAtoms(bool bUnsortedAllowed
, ExceptionInitializer
* errors
) const;
378 * Checks whether the element or its children have unresolved index
383 bool requiresIndexGroups() const;
385 * Resolves an unresolved reference to an index group.
387 * \param[in] grps Index groups to use to resolve the reference.
388 * \param[in] natoms Maximum number of atoms the selections can evaluate to
389 * (zero if the topology/atom count is not set yet).
390 * \throws std::bad_alloc if out of memory.
391 * \throws InconsistentInputError if the reference cannot be
394 void resolveIndexGroupReference(gmx_ana_indexgrps_t
* grps
, int natoms
);
396 * Checks that an index group has valid atom indices.
398 * \param[in] natoms Maximum number of atoms the selections can evaluate to.
399 * \throws std::bad_alloc if out of memory.
400 * \throws InconsistentInputError if there are invalid atom indices.
402 void checkIndexGroup(int natoms
);
404 //! Type of the element.
407 * Value storage of the element.
409 * This field contains the evaluated value of the element, as well as
410 * the output value type.
412 gmx_ana_selvalue_t v
;
414 * Evaluation function for the element.
416 * Can be either NULL (if the expression is a constant and does not
417 * require evaluation) or point to one of the functions defined in
420 sel_evalfunc evaluate
;
422 * Information flags about the element.
424 * Allowed flags are listed here:
425 * \ref selelem_flags "flags for gmx::SelectionTreeElement".
428 //! Data required by the evaluation function.
430 /*! \brief Index group data for several element types.
432 * - \ref SEL_CONST : if the value type is \ref GROUP_VALUE,
433 * this field holds the unprocessed group value.
434 * - \ref SEL_ROOT : holds the group value for which the
435 * selection subtree should be evaluated.
436 * - \ref SEL_SUBEXPR : holds the group for which the subexpression
437 * has been evaluated.
439 gmx_ana_index_t cgrp
;
440 //! Data for \ref SEL_EXPRESSION and \ref SEL_MODIFIER elements.
443 //! Pointer the method used in this expression.
444 struct gmx_ana_selmethod_t
* method
;
445 //! Pointer to the data allocated by the method's \p init_data (see sel_datafunc()).
447 //! Pointer to the position data passed to the method.
448 struct gmx_ana_pos_t
* pos
;
449 //! Pointer to the evaluation data for \p pos.
450 struct gmx_ana_poscalc_t
* pc
;
452 //! Operation type for \ref SEL_BOOLEAN elements.
454 //! Operation type for \ref SEL_ARITHMETIC elements.
459 //! String representation.
462 //! Associated selection parameter for \ref SEL_SUBEXPRREF elements.
463 struct gmx_ana_selparam_t
* param
;
464 //! The string/number used to reference the group.
467 //! Name of the referenced external group.
469 //! If \a name is NULL, the index number of the referenced group.
473 //! Memory pool to use for values, or NULL if standard memory handling.
474 struct gmx_sel_mempool_t
* mempool
;
475 //! Internal data for the selection compiler.
476 t_compiler_data
* cdata
;
478 /*! \brief The first child element.
480 * Other children can be accessed through the \p next field of \p child.
482 SelectionTreeElementPointer child
;
483 //! The next sibling element.
484 SelectionTreeElementPointer next
;
488 * Name of the element.
490 * This field is only used for diagnostic purposes.
494 * Location of the element in the selection text.
496 * This field is only used for diagnostic purposes (including error
499 SelectionLocation location_
;
501 GMX_DISALLOW_COPY_AND_ASSIGN(SelectionTreeElement
);
506 /********************************************************************/
507 /*! \name Selection expression functions
512 /** Writes out a human-readable name for an evaluation function. */
513 void _gmx_sel_print_evalfunc_name(FILE* fp
, gmx::sel_evalfunc evalfunc
);
515 /** Sets the value type of a gmx::SelectionTreeElement. */
516 void _gmx_selelem_set_vtype(const gmx::SelectionTreeElementPointer
& sel
, e_selvalue_t vtype
);
519 * Frees the memory allocated for a selection method parameter.
521 * \param[in] param Parameter to free.
523 void _gmx_selelem_free_param(struct gmx_ana_selparam_t
* param
);
525 * Frees the memory allocated for a selection method.
527 * \param[in] method Method to free.
528 * \param[in] mdata Method data to free.
530 void _gmx_selelem_free_method(struct gmx_ana_selmethod_t
* method
, void* mdata
);
532 /** Prints a human-readable version of a selection element subtree. */
533 void _gmx_selelem_print_tree(FILE* fp
, const gmx::SelectionTreeElement
& sel
, bool bValues
, int level
);
535 /** Prints a human-readable version of the internal compiler data structure. */
536 void _gmx_selelem_print_compiler_info(FILE* fp
, const gmx::SelectionTreeElement
& sel
, int level
);
538 /* In sm_insolidangle.c */
539 /** Returns true if the covered fraction of the selection can be calculated. */
540 bool _gmx_selelem_can_estimate_cover(const gmx::SelectionTreeElement
& sel
);
541 /** Returns the covered fraction of the selection for the current frame. */
542 real
_gmx_selelem_estimate_coverfrac(const gmx::SelectionTreeElement
& sel
);