3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
9 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
10 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
11 * Copyright (c) 2001-2009, The GROMACS development team,
12 * check out http://www.gromacs.org for more information.
14 * This program is free software; you can redistribute it and/or
15 * modify it under the terms of the GNU General Public License
16 * as published by the Free Software Foundation; either version 2
17 * of the License, or (at your option) any later version.
19 * If you want to redistribute modifications, please consider that
20 * scientific software is very special. Version control is crucial -
21 * bugs must be traceable. We will be happy to consider code for
22 * inclusion in the official distribution, but derived work must not
23 * be called official GROMACS. Details are found in the README & COPYING
24 * files - if they are missing, get the official version at www.gromacs.org.
26 * To help us fund GROMACS development, we humbly ask that you cite
27 * the papers on the package - you can find them in the top README file.
29 * For more info, check our website at http://www.gromacs.org
32 * \brief Definition of \c t_selelem and related things.
34 * The selection element trees constructed by the parser and the compiler
35 * are described on the respective pages:
36 * \ref selparser and \ref selcompiler.
38 * This is an implementation header: there should be no need to use it outside
41 #ifndef SELECTION_ELEMENT_H
42 #define SELECTION_ELEMENT_H
44 #include <types/simple.h>
46 #include <indexutil.h>
49 struct gmx_ana_poscalc_t
;
50 struct gmx_ana_selparam_t
;
51 struct gmx_ana_selmethod_t
;
53 struct gmx_sel_evaluate_t
;
54 struct gmx_sel_mempool_t
;
57 /********************************************************************/
58 /*! \name Enumerations for expression types
59 ********************************************************************/
62 /** Defines the type of a \c t_selelem object. */
65 /** Constant-valued expression. */
67 /** Method expression that requires evaluation. */
69 /** Boolean expression. */
71 /** Arithmetic expression. */
73 /** Root node of the evaluation tree. */
75 /** Subexpression that may be referenced several times. */
77 /** Reference to a subexpression. */
79 /** Post-processing of selection value. */
83 /** Defines the boolean operation of \c t_selelem objects with type \ref SEL_BOOLEAN. */
89 BOOL_XOR
/**< Xor (not implemented). */
92 /** Defines the arithmetic operation of \c t_selelem objects with type \ref SEL_ARITHMETIC. */
96 ARITH_MINUS
, /**< - */
97 ARITH_NEG
, /**< Unary - */
100 ARITH_EXP
, /**< ^ (to power) */
103 /** Returns a string representation of the type of a \c t_selelem. */
105 _gmx_selelem_type_str(struct t_selelem
*sel
);
106 /** Returns a string representation of the boolean type of a \ref SEL_BOOLEAN \c t_selelem. */
108 _gmx_selelem_boolean_type_str(struct t_selelem
*sel
);
109 /** Returns a string representation of the type of a \c gmx_ana_selvalue_t. */
111 _gmx_sel_value_type_str(gmx_ana_selvalue_t
*val
);
116 /********************************************************************/
117 /*! \name Selection expression flags
118 * \anchor selelem_flags
119 ********************************************************************/
122 * Selection value flags are set.
124 * If this flag is set, the flags covered by \ref SEL_VALFLAGMASK
125 * have been set properly for the element.
127 #define SEL_FLAGSSET 1
129 * The element evaluates to a single value.
131 * This flag is always set for \ref GROUP_VALUE elements.
133 #define SEL_SINGLEVAL 2
135 * The element evaluates to one value for each input atom.
137 #define SEL_ATOMVAL 4
139 * The element evaluates to an arbitrary number of values.
141 #define SEL_VARNUMVAL 8
143 * The element (or one of its children) is dynamic.
145 #define SEL_DYNAMIC 16
147 * Mask that covers the flags that describe the number of values.
149 #define SEL_VALTYPEMASK (SEL_SINGLEVAL | SEL_ATOMVAL | SEL_VARNUMVAL)
151 * Mask that covers the flags that describe the value type.
153 #define SEL_VALFLAGMASK (SEL_FLAGSSET | SEL_VALTYPEMASK | SEL_DYNAMIC)
155 * Data has been allocated for the \p v.u union.
157 * If not set, the \p v.u.ptr points to data allocated externally.
158 * This is the case if the value of the element is used as a parameter
159 * for a selection method or if the element evaluates the final value of
162 * Even if the flag is set, \p v.u.ptr can be NULL during initialization.
165 * This flag overlaps with the function of \p v.nalloc field, and could
166 * probably be removed, making memory management simpler. Currently, the
167 * \p v.nalloc field is not kept up-to-date in all cases when this flag
168 * is changed and is used in places where this flag is not, so this would
169 * require a careful investigation of the selection code.
171 #define SEL_ALLOCVAL (1<<8)
173 * Data has been allocated for the group/position structure.
175 * If not set, the memory allocated for fields in \p v.u.g or \p v.u.p is
176 * managed externally.
178 * This field has no effect if the value type is not \ref GROUP_VALUE or
179 * \ref POS_VALUE, but should not be set.
181 #define SEL_ALLOCDATA (1<<9)
183 * \p method->init_frame should be called for the frame.
185 #define SEL_INITFRAME (1<<10)
187 * Parameter has been evaluated for the current frame.
189 * This flag is set for children of \ref SEL_EXPRESSION elements (which
190 * describe method parameters) after the element has been evaluated for the
192 * It is not set for \ref SEL_ATOMVAL elements, because they may need to
193 * be evaluated multiple times.
195 #define SEL_EVALFRAME (1<<11)
197 * \p method->init has been called.
199 #define SEL_METHODINIT (1<<12)
201 * \p method->outinit has been called.
203 * This flag is also used for \ref SEL_SUBEXPRREF elements.
205 #define SEL_OUTINIT (1<<13)
209 /********************************************************************/
210 /*! \name Selection expression data structures and functions
211 ********************************************************************/
217 * Function pointer for evaluating a \c t_selelem.
219 typedef int (*sel_evalfunc
)(struct gmx_sel_evaluate_t
*data
,
220 struct t_selelem
*sel
, gmx_ana_index_t
*g
);
223 * Represents an element of a selection expression.
225 typedef struct t_selelem
227 /*! \brief Name of the element.
229 * This field is only used for informative purposes.
230 * It is always either NULL or a pointer to a string.
231 * Memory is never allocated for it directly.
234 /** Type of the element. */
237 * Value storage of the element.
239 * This field contains the evaluated value of the element, as well as
240 * the output value type.
242 gmx_ana_selvalue_t v
;
244 * Evaluation function for the element.
246 * Can be either NULL (if the expression is a constant and does not require
247 * evaluation) or point to one of the functions defined in evaluate.h.
249 sel_evalfunc evaluate
;
251 * Information flags about the element.
253 * Allowed flags are listed here:
254 * \ref selelem_flags "flags for \c t_selelem".
257 /** Data required by the evaluation function. */
259 /*! \brief Index group data for several element types.
261 * - \ref SEL_CONST : if the value type is \ref GROUP_VALUE,
262 * this field holds the unprocessed group value.
263 * - \ref SEL_ROOT : holds the group value for which the
264 * selection subtree should be evaluated.
265 * - \ref SEL_SUBEXPR : holds the group for which the subexpression
266 * has been evaluated.
268 gmx_ana_index_t cgrp
;
269 /** Data for \ref SEL_EXPRESSION and \ref SEL_MODIFIER elements. */
271 /** Pointer the the method used in this expression. */
272 struct gmx_ana_selmethod_t
*method
;
273 /** Pointer to the data allocated by the method's \p init_data (see sel_datafunc()). */
275 /** Pointer to the position data passed to the method. */
276 struct gmx_ana_pos_t
*pos
;
277 /** Pointer to the evaluation data for \p pos. */
278 struct gmx_ana_poscalc_t
*pc
;
280 /** Operation type for \ref SEL_BOOLEAN elements. */
282 /** Operation type for \ref SEL_ARITHMETIC elements. */
284 /** Operation type. */
286 /** String representation. */
289 /** Associated selection parameter for \ref SEL_SUBEXPRREF elements. */
290 struct gmx_ana_selparam_t
*param
;
292 /** Memory pool to use for values, or NULL if standard memory handling. */
293 struct gmx_sel_mempool_t
*mempool
;
294 /** Internal data for the selection compiler. */
295 struct t_compiler_data
*cdata
;
297 /*! \brief The first child element.
299 * Other children can be accessed through the \p next field of \p child.
301 struct t_selelem
*child
;
302 /** The next sibling element. */
303 struct t_selelem
*next
;
304 /*! \brief Number of references to this element.
306 * Should be larger than one only for \ref SEL_SUBEXPR elements.
311 /** Allocates memory and performs some common initialization for a \c t_selelem. */
313 _gmx_selelem_create(e_selelem_t type
);
314 /** Sets the value type of a \c t_selelem. */
316 _gmx_selelem_set_vtype(t_selelem
*sel
, e_selvalue_t vtype
);
317 /** Reserves memory for value of a \c t_selelem from a memory pool. */
319 _gmx_selelem_mempool_reserve(t_selelem
*sel
, int count
);
320 /** Releases memory pool used for value of a \c t_selelem. */
322 _gmx_selelem_mempool_release(t_selelem
*sel
);
323 /** Frees the memory allocated for a \c t_selelem structure and all its children. */
325 _gmx_selelem_free(t_selelem
*sel
);
326 /** Frees the memory allocated for a \c t_selelem structure, all its children, and also all structures referenced through t_selelem::next fields. */
328 _gmx_selelem_free_chain(t_selelem
*first
);
330 /** Frees the memory allocated for the \c t_selelem::d union. */
332 _gmx_selelem_free_values(t_selelem
*sel
);
333 /** Frees the memory allocated for a selection method. */
335 _gmx_selelem_free_method(struct gmx_ana_selmethod_t
*method
, void *mdata
);
336 /** Frees the memory allocated for the \c t_selelem::u field. */
338 _gmx_selelem_free_exprdata(t_selelem
*sel
);
340 /** Frees the memory allocated for the selection compiler. */
342 _gmx_selelem_free_compiler_data(t_selelem
*sel
);
344 /** Prints a human-readable version of a selection element subtree. */
346 _gmx_selelem_print_tree(FILE *fp
, t_selelem
*root
, bool bValues
, int level
);
348 /** Returns TRUE if the selection element subtree requires topology information for evaluation. */
350 _gmx_selelem_requires_top(t_selelem
*root
);
352 /* In sm_insolidangle.c */
353 /** Returns TRUE if the covered fraction of the selection can be calculated. */
355 _gmx_selelem_can_estimate_cover(t_selelem
*sel
);
356 /** Returns the covered fraction of the selection for the current frame. */
358 _gmx_selelem_estimate_coverfrac(t_selelem
*sel
);