Merge branch release-4-6
[alexxy/gromacs.git] / src / gromacs / selection / poscalc.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2009,2010,2011,2012, by the GROMACS development team, led by
5  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
6  * others, as listed in the AUTHORS file in the top-level source
7  * 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  * API for structured and optimized calculation of positions.
38  *
39  * This header declares an API for calculating positions in an automated way,
40  * for internal use by the selection engine.  This is useful in particular with
41  * dynamic selections, because the same COM/COG positions may be needed in
42  * several contexts.  The API makes it possible to optimize the evaluation such
43  * that any heavy calculation is only done once, and the results just copied if
44  * needed more than once.  The functions also provide a convenient interface
45  * for keeping the whole \c gmx_ana_pos_t structure up-to-date.
46  *
47  * The API is documented in more detail in gmx::PositionCalculationCollection.
48  *
49  * \author Teemu Murtola <teemu.murtola@gmail.com>
50  * \ingroup module_selection
51  */
52 #ifndef GMX_SELECTION_POSCALC_H
53 #define GMX_SELECTION_POSCALC_H
54
55 #include "../legacyheaders/typedefs.h"
56
57 #include "../utility/common.h"
58
59 /*! \name Flags for position calculation.
60  * \anchor poscalc_flags
61  */
62 /*@{*/
63 /*! \brief
64  * Use mass weighting.
65  *
66  * If this flag is set, the positions will be calculated using mass weighting,
67  * i.e., one gets center-of-mass positions.
68  * Without the flag, center-of-geometry positions are calculated.
69  * Does not have any effect if the calculation type is \ref POS_ATOM.
70  */
71 #define POS_MASS        1
72 /*! \brief
73  * Calculate positions for the same atoms in residues/molecules.
74  *
75  * If this flag is set, the positions are always calculated using the same
76  * atoms for each residue/molecule, even if the evaluation group contains only
77  * some of the atoms for some frames.
78  * The group passed to gmx_ana_poscalc_set_maxindex() is used to determine
79  * the atoms to use for the calculation.
80  *
81  * Has no effect unless \ref POS_DYNAMIC is set or if the calculation type
82  * is not \ref POS_RES of \ref POS_MOL.
83  */
84 #define POS_COMPLMAX    2
85 /*! \brief
86  * Calculate positions for whole residues/molecules.
87  *
88  * If this flag is set, the positions will be calculated for whole
89  * residues/molecules, even if the group contains only some of the atoms in
90  * the residue/molecule.
91  *
92  * Has no effect unless the calculation type is \ref POS_RES or \ref POS_MOL.
93  */
94 #define POS_COMPLWHOLE  4
95 /*! \brief
96  * Enable handling of changing calculation groups.
97  *
98  * Can be used for static calculations as well, but implies a small
99  * performance penalty.
100  */
101 #define POS_DYNAMIC     16
102 /*! \brief
103  * Update \c gmx_ana_pos_t::m dynamically for an otherwise static
104  * calculation.
105  *
106  * Has effect only if \ref POS_DYNAMIC is not set.
107  */
108 #define POS_MASKONLY    32
109 /*! \brief
110  * Calculate velocities of the positions.
111  */
112 #define POS_VELOCITIES  64
113 /*! \brief
114  * Calculate forces on the positions.
115  */
116 #define POS_FORCES      128
117 /*@}*/
118
119 /** Specifies the type of positions to be calculated. */
120 typedef enum
121 {
122     POS_ATOM,    /**< Copy atomic coordinates. */
123     POS_RES,     /**< Calculate center for each residue. */
124     POS_MOL,     /**< Calculate center for each molecule. */
125     POS_ALL,     /**< Calculate center for the whole group. */
126     POS_ALL_PBC  /**< Calculate center for the whole group with PBC. */
127 } e_poscalc_t;
128
129 /** Data structure for position calculation. */
130 typedef struct gmx_ana_poscalc_t gmx_ana_poscalc_t;
131
132 struct gmx_ana_index_t;
133 struct gmx_ana_pos_t;
134
135 namespace gmx
136 {
137
138 /*! \internal \brief
139  * Collection of \c gmx_ana_poscalc_t structures for the same topology.
140  *
141  * Calculations within one collection share the same topology, and they are
142  * optimized.  Calculations in different collections do not interact.
143  * The topology for a collection can be set with setTopology().
144  * This needs to be done before calling gmx_ana_poscalc_set_maxindex() for
145  * any calculation in the collection, unless that calculation does not
146  * require topology information.
147  *
148  * A new calculation is created with createCalculation().
149  * If flags need to be adjusted later, gmx_ana_poscalc_set_flags() can be
150  * used.
151  * After the flags are final, the largest possible index group for which the
152  * positions are needed has to be set with gmx_ana_poscalc_set_maxindex().
153  * setTopology() should have been called before this function is called.
154  * After the above calls, gmx_ana_poscalc_init_pos() can be used to initialize
155  * output to a \c gmx_ana_pos_t structure.  Several different structures can be
156  * initialized for the same calculation; the only requirement is that the
157  * structure passed later to gmx_ana_poscalc_update() has been initialized
158  * properly.
159  * The memory allocated for a calculation can be freed with
160  * gmx_ana_poscalc_free().
161  *
162  * The position evaluation is simple: initFrame() should be
163  * called once for each frame, and gmx_ana_poscalc_update() can then be called
164  * for each calculation that is needed for that frame.
165  *
166  * It is also possible to initialize the calculations based on a type provided
167  * as a string.
168  * The possible strings are listed in \ref typeEnumValues, and the string can
169  * be converted to the parameters for createCalculation() using typeFromEnum().
170  * createCalculationFromEnum() is also provided for convenience.
171  *
172  * \ingroup module_selection
173  */
174 class PositionCalculationCollection
175 {
176     public:
177         /*! \brief
178          * Array of strings acceptable for position calculation type enum.
179          *
180          * This array contains the acceptable values for typeFromEnum() and
181          * createCalculationFromEnum().
182          * The array contains a NULL pointer after the last item to indicate
183          * the end of the list.
184          */
185         static const char * const typeEnumValues[];
186
187         /*! \brief
188          * Converts a string to parameters for createCalculationFromEnum().
189          *
190          * \param[in]     post  String (typically an enum argument).
191          *     Allowed values: 'atom', 'res_com', 'res_cog', 'mol_com', 'mol_cog',
192          *     or one of the last four prepended by 'whole_', 'part_', or 'dyn_'.
193          * \param[out]    type  \c e_poscalc_t corresponding to \p post.
194          * \param[in,out] flags Flags corresponding to \p post.
195          *     On input, the flags should contain the default flags.
196          *     On exit, the flags \ref POS_MASS, \ref POS_COMPLMAX and
197          *     \ref POS_COMPLWHOLE have been set according to \p post
198          *     (the completion flags are left at the default values if no
199          *     completion prefix is given).
200          * \throws  InternalError  if post is not recognized.
201          *
202          * \attention
203          * Checking is not complete, and other values than those listed above
204          * may be accepted for \p post, but the results are undefined.
205          *
206          * \see typeEnumValues
207          */
208         static void typeFromEnum(const char *post, e_poscalc_t *type, int *flags);
209
210         /*! \brief
211          * Creates a new position calculation collection object.
212          *
213          * \throws  std::bad_alloc if out of memory.
214          */
215         PositionCalculationCollection();
216         /*! \brief
217          * Destroys a position calculation collection and its calculations.
218          *
219          * Any calculations in the collection are also freed, even if
220          * references to them are left.
221          */
222         ~PositionCalculationCollection();
223
224         /*! \brief
225          * Sets the topology used for the calculations.
226          *
227          * \param[in]     top   Topology data structure.
228          *
229          * This function should be called to set the topology before using
230          * gmx_ana_poscalc_set_maxindex() for any calculation that requires
231          * topology information.
232          *
233          * Does not throw.
234          */
235         void setTopology(t_topology *top);
236         /*! \brief
237          * Prints information about calculations.
238          *
239          * \param[in] fp    File handle to receive the output.
240          *
241          * The output is very technical, making this function mainly useful for
242          * debugging purposes.
243          *
244          * Does not throw.
245          */
246         void printTree(FILE *fp) const;
247
248         /*! \brief
249          * Creates a new position calculation.
250          *
251          * \param[in]  type  Type of calculation.
252          * \param[in]  flags Flags for setting calculation options
253          *   (see \ref poscalc_flags "documentation of the flags").
254          *
255          * Does not throw currently, but may throw std::bad_alloc in the
256          * future.
257          */
258         gmx_ana_poscalc_t *createCalculation(e_poscalc_t type, int flags);
259         /*! \brief
260          * Creates a new position calculation based on an enum value.
261          *
262          * \param[in]  post  One of the strings acceptable for
263          *      typeFromEnum().
264          * \param[in]  flags Flags for setting calculation options
265          *      (see \ref poscalc_flags "documentation of the flags").
266          * \throws     InternalError  if post is not recognized.
267          *
268          * This is a convenience wrapper for createCalculation().
269          * \p flags sets the default calculation options if not overridden by
270          * \p post; see typeFromEnum().
271          *
272          * May also throw std::bad_alloc in the future.
273          *
274          * \see createCalculation(), typeFromEnum()
275          */
276         gmx_ana_poscalc_t *createCalculationFromEnum(const char *post, int flags);
277
278         /*! \brief
279          * Initializes evaluation for a position calculation collection.
280          *
281          * This function does some final initialization of the data structures
282          * in the collection to prepare them for evaluation.
283          * After this function has been called, it is no longer possible to add
284          * new calculations to the collection.
285          *
286          * Multiple calls to the function are ignored.
287          *
288          * Does not throw currently, but may throw std::bad_alloc in the
289          * future.
290          */
291         void initEvaluation();
292         /*! \brief
293          * Initializes a position calculation collection for a new frame.
294          *
295          * Clears the evaluation flag for all calculations.
296          * Should be called for each frame before calling
297          * gmx_ana_poscalc_update().
298          *
299          * This function calls initEvaluation() automatically if it has not
300          * been called earlier.
301          *
302          * Does not throw, but may throw if initEvaluation() is changed to
303          * throw.
304          */
305         void initFrame();
306
307     private:
308         class Impl;
309
310         PrivateImplPointer<Impl> impl_;
311
312         /*! \brief
313          * Needed to access the implementation class from the C code.
314          */
315         friend struct ::gmx_ana_poscalc_t;
316 };
317
318 } // namespace gmx
319
320 /** Sets the flags for position calculation. */
321 void
322 gmx_ana_poscalc_set_flags(gmx_ana_poscalc_t *pc, int flags);
323 /** Sets the maximum possible input index group for position calculation. */
324 void
325 gmx_ana_poscalc_set_maxindex(gmx_ana_poscalc_t *pc, struct gmx_ana_index_t *g);
326 /** Initializes positions for position calculation output. */
327 void
328 gmx_ana_poscalc_init_pos(gmx_ana_poscalc_t *pc, struct gmx_ana_pos_t *p);
329 /** Frees the memory allocated for position calculation. */
330 void
331 gmx_ana_poscalc_free(gmx_ana_poscalc_t *pc);
332 /** Returns true if the position calculation requires topology information. */
333 bool
334 gmx_ana_poscalc_requires_top(gmx_ana_poscalc_t *pc);
335
336 /** Updates a single COM/COG structure for a frame. */
337 void
338 gmx_ana_poscalc_update(gmx_ana_poscalc_t *pc,
339                        struct gmx_ana_pos_t *p, struct gmx_ana_index_t *g,
340                        t_trxframe *fr, t_pbc *pbc);
341
342 #endif