1a1e86bf611d7644c65de9dcf8748dfb27d0f97e
[alexxy/gromacs.git] / src / gromacs / selection / selelem.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2009,2010,2011,2012,2013,2014, 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  * Declares gmx::SelectionTreeElement and related things.
38  *
39  * The selection element trees constructed by the parser and the compiler
40  * are described on the respective pages:
41  * \ref page_module_selection_parser and \ref page_module_selection_compiler.
42  *
43  * This is an implementation header: there should be no need to use it outside
44  * this directory.
45  *
46  * \author Teemu Murtola <teemu.murtola@gmail.com>
47  * \ingroup module_selection
48  */
49 #ifndef GMX_SELECTION_SELELEM_H
50 #define GMX_SELECTION_SELELEM_H
51
52 #include <string>
53
54 #include <boost/shared_ptr.hpp>
55
56 #include "gromacs/utility/common.h"
57 #include "gromacs/utility/real.h"
58
59 #include "indexutil.h"
60 #include "selvalue.h"
61
62 struct gmx_ana_poscalc_t;
63 struct gmx_ana_selparam_t;
64 struct gmx_ana_selmethod_t;
65
66 struct gmx_sel_evaluate_t;
67 struct gmx_sel_mempool_t;
68
69 struct t_compiler_data;
70
71 namespace gmx
72 {
73 class SelectionTreeElement;
74
75 //! Smart pointer type for selection tree element pointers.
76 typedef boost::shared_ptr<SelectionTreeElement> SelectionTreeElementPointer;
77 } // namespace gmx
78
79 /********************************************************************/
80 /*! \name Enumerations for expression types
81  ********************************************************************/
82 //!\{
83
84 /** Defines the type of a gmx::SelectionTreeElement object. */
85 typedef enum
86 {
87     /** Constant-valued expression. */
88     SEL_CONST,
89     /** Method expression that requires evaluation. */
90     SEL_EXPRESSION,
91     /** Boolean expression. */
92     SEL_BOOLEAN,
93     /** Arithmetic expression. */
94     SEL_ARITHMETIC,
95     /** Root node of the evaluation tree. */
96     SEL_ROOT,
97     /** Subexpression that may be referenced several times. */
98     SEL_SUBEXPR,
99     /** Reference to a subexpression. */
100     SEL_SUBEXPRREF,
101     /** Unresolved reference to an external group. */
102     SEL_GROUPREF,
103     /** Post-processing of selection value. */
104     SEL_MODIFIER
105 } e_selelem_t;
106
107 /** Defines the boolean operation of gmx::SelectionTreeElement objects with type \ref SEL_BOOLEAN. */
108 typedef enum
109 {
110     BOOL_NOT,           /**< Not */
111     BOOL_AND,           /**< And */
112     BOOL_OR,            /**< Or */
113     BOOL_XOR            /**< Xor (not implemented). */
114 } e_boolean_t;
115
116 /** Defines the arithmetic operation of gmx::SelectionTreeElement objects with type \ref SEL_ARITHMETIC. */
117 typedef enum
118 {
119     ARITH_PLUS,         /**< Addition (`+`) */
120     ARITH_MINUS,        /**< Subtraction (`-`) */
121     ARITH_NEG,          /**< Unary `-` */
122     ARITH_MULT,         /**< Multiplication (`*`) */
123     ARITH_DIV,          /**< Division (`/`) */
124     ARITH_EXP           /**< Power (`^`) */
125 } e_arithmetic_t;
126
127 /** Returns a string representation of the type of a gmx::SelectionTreeElement. */
128 extern const char *
129 _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 *
132 _gmx_selelem_boolean_type_str(const gmx::SelectionTreeElement &sel);
133 /** Returns a string representation of the type of a \c gmx_ana_selvalue_t. */
134 extern const char *
135 _gmx_sel_value_type_str(const gmx_ana_selvalue_t *val);
136
137 //!\}
138
139
140 /********************************************************************/
141 /*! \name Selection expression flags
142  * \anchor selelem_flags
143  ********************************************************************/
144 //!\{
145 /*! \brief
146  * Selection value flags are set.
147  *
148  * If this flag is set, the flags covered by \ref SEL_VALFLAGMASK
149  * have been set properly for the element.
150  */
151 #define SEL_FLAGSSET    1
152 /*! \brief
153  * The element evaluates to a single value.
154  *
155  * This flag is always set for \ref GROUP_VALUE elements.
156  */
157 #define SEL_SINGLEVAL   2
158 /*! \brief
159  * The element evaluates to one value for each input atom.
160  */
161 #define SEL_ATOMVAL     4
162 /*! \brief
163  * The element evaluates to an arbitrary number of values.
164  */
165 #define SEL_VARNUMVAL   8
166 /*! \brief
167  * The element (or one of its children) is dynamic.
168  */
169 #define SEL_DYNAMIC     16
170 /*! \brief
171  * The element may contain atom indices in an unsorted order.
172  */
173 #define SEL_UNSORTED    32
174 /*! \brief
175  * Mask that covers the flags that describe the number of values.
176  */
177 #define SEL_VALTYPEMASK (SEL_SINGLEVAL | SEL_ATOMVAL | SEL_VARNUMVAL)
178 /*! \brief
179  * Mask that covers the flags that describe the value type.
180  */
181 #define SEL_VALFLAGMASK (SEL_FLAGSSET | SEL_VALTYPEMASK | SEL_DYNAMIC)
182 /*! \brief
183  * Data has been allocated for the \p v.u union.
184  *
185  * If not set, the \p v.u.ptr points to data allocated externally.
186  * This is the case if the value of the element is used as a parameter
187  * for a selection method or if the element evaluates the final value of
188  * a selection.
189  *
190  * Even if the flag is set, \p v.u.ptr can be NULL during initialization.
191  *
192  * \todo
193  * This flag overlaps with the function of \p v.nalloc field, and could
194  * probably be removed, making memory management simpler. Currently, the
195  * \p v.nalloc field is not kept up-to-date in all cases when this flag
196  * is changed and is used in places where this flag is not, so this would
197  * require a careful investigation of the selection code.
198  */
199 #define SEL_ALLOCVAL    (1<<8)
200 /*! \brief
201  * Data has been allocated for the group/position structure.
202  *
203  * If not set, the memory allocated for fields in \p v.u.g or \p v.u.p is
204  * managed externally.
205  *
206  * This field has no effect if the value type is not \ref GROUP_VALUE or
207  * \ref POS_VALUE, but should not be set.
208  */
209 #define SEL_ALLOCDATA   (1<<9)
210 /*! \brief
211  * \p method->init_frame should be called for the frame.
212  */
213 #define SEL_INITFRAME   (1<<10)
214 /*! \brief
215  * Parameter has been evaluated for the current frame.
216  *
217  * This flag is set for children of \ref SEL_EXPRESSION elements (which
218  * describe method parameters) after the element has been evaluated for the
219  * current frame.
220  * It is not set for \ref SEL_ATOMVAL elements, because they may need to
221  * be evaluated multiple times.
222  */
223 #define SEL_EVALFRAME   (1<<11)
224 /*! \brief
225  * \p method->init has been called.
226  */
227 #define SEL_METHODINIT  (1<<12)
228 /*! \brief
229  * \p method->outinit has been called.
230  *
231  * This flag is also used for \ref SEL_SUBEXPRREF elements.
232  */
233 #define SEL_OUTINIT     (1<<13)
234 //!\}
235
236
237 namespace gmx
238 {
239
240 class ExceptionInitializer;
241
242 /*! \brief
243  * Function pointer for evaluating a gmx::SelectionTreeElement.
244  */
245 typedef void (*sel_evalfunc)(struct gmx_sel_evaluate_t         *data,
246                              const SelectionTreeElementPointer &sel,
247                              gmx_ana_index_t                   *g);
248
249 /*! \internal \brief
250  * Represents an element of a selection expression.
251  */
252 class SelectionTreeElement
253 {
254     public:
255         /*! \brief
256          * Allocates memory and performs common initialization.
257          *
258          * \param[in] type Type of selection element to create.
259          *
260          * \a type is set to \p type,
261          * \a v::type is set to \ref GROUP_VALUE for boolean and comparison
262          * expressions and \ref NO_VALUE for others, and
263          * \ref SEL_ALLOCVAL is set for non-root elements (\ref SEL_ALLOCDATA
264          * is also set for \ref SEL_BOOLEAN elements).
265          * All the pointers are set to NULL.
266          */
267         explicit SelectionTreeElement(e_selelem_t type);
268         ~SelectionTreeElement();
269
270         //! Frees the memory allocated for the \a v union.
271         void freeValues();
272         //! Frees the memory allocated for the \a u union.
273         void freeExpressionData();
274         /* In compiler.cpp */
275         /*! \brief
276          * Frees the memory allocated for the selection compiler.
277          *
278          * This function only frees the data for the given selection, not its
279          * children.  It is safe to call the function when compiler data has
280          * not been allocated or has already been freed; in such a case,
281          * nothing is done.
282          */
283         void freeCompilerData();
284
285         /*! \brief
286          * Reserves memory for value from a memory pool.
287          *
288          * \param[in]     count Number of values to reserve memory for.
289          *
290          * Reserves memory for the values of this element from the \a mempool
291          * memory pool.
292          * If no memory pool is set, nothing is done.
293          */
294         void mempoolReserve(int count);
295         /*! \brief
296          * Releases memory pool used for value.
297          *
298          * Releases the memory allocated for the values of this element from the
299          * \a mempool memory pool.
300          * If no memory pool is set, nothing is done.
301          */
302         void mempoolRelease();
303
304         //! Returns the name of the element.
305         const std::string &name() const { return name_; }
306         /*! \brief
307          * Sets the name of the element.
308          *
309          * \param[in] name  Name to set (can be NULL).
310          * \throws    std::bad_alloc if out of memory.
311          */
312         void setName(const char *name) { name_ = (name != NULL ? name : ""); }
313         //! \copydoc setName(const char *)
314         void setName(const std::string &name) { name_ = name; }
315         /*! \brief
316          * Sets the name of a root element if it is missing.
317          *
318          * \param[in] selectionText  Full selection text to use as a fallback.
319          * \throws    std::bad_alloc if out of memory.
320          *
321          * If index groups have not yet been set and the selection is a result
322          * of a group reference, the name may still be empty after this call.
323          *
324          * Strong exception safety guarantee.
325          */
326         void fillNameIfMissing(const char *selectionText);
327
328         /*! \brief
329          * Checks that this element and its children do not contain unsupported
330          * elements with unsorted atoms.
331          *
332          * \param[in] bUnsortedAllowed Whether this element's parents allow it
333          *     to have unsorted atoms.
334          * \param     errors           Object for reporting any error messages.
335          * \throws    std::bad_alloc if out of memory.
336          *
337          * Errors are reported as nested exceptions in \p errors.
338          */
339         void checkUnsortedAtoms(bool                  bUnsortedAllowed,
340                                 ExceptionInitializer *errors) const;
341         /*! \brief
342          * Resolves an unresolved reference to an index group.
343          *
344          * \param[in] grps   Index groups to use to resolve the reference.
345          * \param[in] natoms Maximum number of atoms the selections can evaluate to
346          *     (zero if the topology/atom count is not set yet).
347          * \throws    std::bad_alloc if out of memory.
348          * \throws    InconsistentInputError if the reference cannot be
349          *     resolved.
350          */
351         void resolveIndexGroupReference(gmx_ana_indexgrps_t *grps, int natoms);
352         /*! \brief
353          * Checks that an index group has valid atom indices.
354          *
355          * \param[in] natoms Maximum number of atoms the selections can evaluate to.
356          * \throws    std::bad_alloc if out of memory.
357          * \throws    InconsistentInputError if there are invalid atom indices.
358          */
359         void checkIndexGroup(int natoms);
360
361         //! Type of the element.
362         e_selelem_t                         type;
363         /*! \brief
364          * Value storage of the element.
365          *
366          * This field contains the evaluated value of the element, as well as
367          * the output value type.
368          */
369         gmx_ana_selvalue_t                  v;
370         /*! \brief
371          * Evaluation function for the element.
372          *
373          * Can be either NULL (if the expression is a constant and does not
374          * require evaluation) or point to one of the functions defined in
375          * evaluate.h.
376          */
377         sel_evalfunc                        evaluate;
378         /*! \brief
379          * Information flags about the element.
380          *
381          * Allowed flags are listed here:
382          * \ref selelem_flags "flags for gmx::SelectionTreeElement".
383          */
384         int                                 flags;
385         //! Data required by the evaluation function.
386         union {
387             /*! \brief Index group data for several element types.
388              *
389              *  - \ref SEL_CONST : if the value type is \ref GROUP_VALUE,
390              *    this field holds the unprocessed group value.
391              *  - \ref SEL_ROOT : holds the group value for which the
392              *    selection subtree should be evaluated.
393              *  - \ref SEL_SUBEXPR : holds the group for which the subexpression
394              *    has been evaluated.
395              */
396             gmx_ana_index_t                 cgrp;
397             //! Data for \ref SEL_EXPRESSION and \ref SEL_MODIFIER elements.
398             struct {
399                 //! Pointer the the method used in this expression.
400                 struct gmx_ana_selmethod_t *method;
401                 //! Pointer to the data allocated by the method's \p init_data (see sel_datafunc()).
402                 void                       *mdata;
403                 //! Pointer to the position data passed to the method.
404                 struct gmx_ana_pos_t       *pos;
405                 //! Pointer to the evaluation data for \p pos.
406                 struct gmx_ana_poscalc_t   *pc;
407             }                               expr;
408             //! Operation type for \ref SEL_BOOLEAN elements.
409             e_boolean_t                     boolt;
410             //! Operation type for \ref SEL_ARITHMETIC elements.
411             struct {
412                 //! Operation type.
413                 e_arithmetic_t              type;
414                 //! String representation.
415                 char                       *opstr;
416             }                               arith;
417             //! Associated selection parameter for \ref SEL_SUBEXPRREF elements.
418             struct gmx_ana_selparam_t      *param;
419             //! The string/number used to reference the group.
420             struct {
421                 //! Name of the referenced external group.
422                 char                       *name;
423                 //! If \a name is NULL, the index number of the referenced group.
424                 int                         id;
425             }                               gref;
426         }                                   u;
427         //! Memory pool to use for values, or NULL if standard memory handling.
428         struct gmx_sel_mempool_t           *mempool;
429         //! Internal data for the selection compiler.
430         t_compiler_data                    *cdata;
431
432         /*! \brief The first child element.
433          *
434          * Other children can be accessed through the \p next field of \p child.
435          */
436         SelectionTreeElementPointer         child;
437         //! The next sibling element.
438         SelectionTreeElementPointer         next;
439
440     private:
441         /*! \brief
442          * Name of the element.
443          *
444          * This field is only used for informative purposes.
445          */
446         std::string                         name_;
447
448         GMX_DISALLOW_COPY_AND_ASSIGN(SelectionTreeElement);
449 };
450
451 } // namespace gmx
452
453 /********************************************************************/
454 /*! \name Selection expression functions
455  */
456 //!\{
457
458 /* In evaluate.c */
459 /** Writes out a human-readable name for an evaluation function. */
460 void
461 _gmx_sel_print_evalfunc_name(FILE *fp, gmx::sel_evalfunc evalfunc);
462
463 /** Sets the value type of a gmx::SelectionTreeElement. */
464 void
465 _gmx_selelem_set_vtype(const gmx::SelectionTreeElementPointer &sel,
466                        e_selvalue_t                            vtype);
467
468 /** Frees the memory allocated for a selection method. */
469 void
470 _gmx_selelem_free_method(struct gmx_ana_selmethod_t *method, void *mdata);
471
472 /** Prints a human-readable version of a selection element subtree. */
473 void
474 _gmx_selelem_print_tree(FILE *fp, const gmx::SelectionTreeElement &sel,
475                         bool bValues, int level);
476 /* In compiler.c */
477 /** Prints a human-readable version of the internal compiler data structure. */
478 void
479 _gmx_selelem_print_compiler_info(FILE *fp, const gmx::SelectionTreeElement &sel,
480                                  int level);
481
482 /** Returns true if the selection element subtree requires topology information for evaluation. */
483 bool
484 _gmx_selelem_requires_top(const gmx::SelectionTreeElement &root);
485
486 /* In sm_insolidangle.c */
487 /** Returns true if the covered fraction of the selection can be calculated. */
488 bool
489 _gmx_selelem_can_estimate_cover(const gmx::SelectionTreeElement &sel);
490 /** Returns the covered fraction of the selection for the current frame. */
491 real
492 _gmx_selelem_estimate_coverfrac(const gmx::SelectionTreeElement &sel);
493
494 //!\}
495
496 #endif