Sort all includes in src/gromacs
[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,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  * 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 <cstdio>
56
57 #include "gromacs/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 struct gmx_ana_poscalc_t;
131
132 struct gmx_ana_index_t;
133 struct gmx_ana_pos_t;
134 struct t_pbc;
135 struct t_topology;
136 struct t_trxframe;
137
138 namespace gmx
139 {
140
141 /*! \internal
142  * \brief
143  * Collection of \c gmx_ana_poscalc_t structures for the same topology.
144  *
145  * Calculations within one collection share the same topology, and they are
146  * optimized.  Calculations in different collections do not interact.
147  * The topology for a collection can be set with setTopology().
148  * This needs to be done before calling gmx_ana_poscalc_set_maxindex() for
149  * any calculation in the collection, unless that calculation does not
150  * require topology information.
151  *
152  * A new calculation is created with createCalculation().
153  * If flags need to be adjusted later, gmx_ana_poscalc_set_flags() can be
154  * used.
155  * After the flags are final, the largest possible index group for which the
156  * positions are needed has to be set with gmx_ana_poscalc_set_maxindex().
157  * setTopology() should have been called before this function is called.
158  * After the above calls, gmx_ana_poscalc_init_pos() can be used to initialize
159  * output to a \c gmx_ana_pos_t structure.  Several different structures can be
160  * initialized for the same calculation; the only requirement is that the
161  * structure passed later to gmx_ana_poscalc_update() has been initialized
162  * properly.
163  * The memory allocated for a calculation can be freed with
164  * gmx_ana_poscalc_free().
165  *
166  * The position evaluation is simple: initFrame() should be
167  * called once for each frame, and gmx_ana_poscalc_update() can then be called
168  * for each calculation that is needed for that frame.
169  *
170  * It is also possible to initialize the calculations based on a type provided
171  * as a string.
172  * The possible strings are listed in \ref typeEnumValues, and the string can
173  * be converted to the parameters for createCalculation() using typeFromEnum().
174  * createCalculationFromEnum() is also provided for convenience.
175  *
176  * \ingroup module_selection
177  */
178 class PositionCalculationCollection
179 {
180     public:
181         /*! \brief
182          * Array of strings acceptable for position calculation type enum.
183          *
184          * This array contains the acceptable values for typeFromEnum() and
185          * createCalculationFromEnum().
186          * The array contains a NULL pointer after the last item to indicate
187          * the end of the list.
188          */
189         static const char * const typeEnumValues[];
190
191         /*! \brief
192          * Converts a string to parameters for createCalculationFromEnum().
193          *
194          * \param[in]     post  String (typically an enum argument).
195          *     Allowed values: 'atom', 'res_com', 'res_cog', 'mol_com', 'mol_cog',
196          *     or one of the last four prepended by 'whole_', 'part_', or 'dyn_'.
197          * \param[out]    type  \c e_poscalc_t corresponding to \p post.
198          * \param[in,out] flags Flags corresponding to \p post.
199          *     On input, the flags should contain the default flags.
200          *     On exit, the flags \ref POS_MASS, \ref POS_COMPLMAX and
201          *     \ref POS_COMPLWHOLE have been set according to \p post
202          *     (the completion flags are left at the default values if no
203          *     completion prefix is given).
204          * \throws  InternalError  if post is not recognized.
205          *
206          * \attention
207          * Checking is not complete, and other values than those listed above
208          * may be accepted for \p post, but the results are undefined.
209          *
210          * \see typeEnumValues
211          */
212         static void typeFromEnum(const char *post, e_poscalc_t *type, int *flags);
213
214         /*! \brief
215          * Creates a new position calculation collection object.
216          *
217          * \throws  std::bad_alloc if out of memory.
218          */
219         PositionCalculationCollection();
220         /*! \brief
221          * Destroys a position calculation collection and its calculations.
222          *
223          * Any calculations in the collection are also freed, even if
224          * references to them are left.
225          */
226         ~PositionCalculationCollection();
227
228         /*! \brief
229          * Sets the topology used for the calculations.
230          *
231          * \param[in]     top   Topology data structure.
232          *
233          * This function should be called to set the topology before using
234          * gmx_ana_poscalc_set_maxindex() for any calculation that requires
235          * topology information.
236          *
237          * Does not throw.
238          */
239         void setTopology(t_topology *top);
240         /*! \brief
241          * Prints information about calculations.
242          *
243          * \param[in] fp    File handle to receive the output.
244          *
245          * The output is very technical, making this function mainly useful for
246          * debugging purposes.
247          *
248          * Does not throw.
249          */
250         void printTree(FILE *fp) const;
251
252         /*! \brief
253          * Creates a new position calculation.
254          *
255          * \param[in]  type  Type of calculation.
256          * \param[in]  flags Flags for setting calculation options
257          *   (see \ref poscalc_flags "documentation of the flags").
258          *
259          * Does not throw currently, but may throw std::bad_alloc in the
260          * future.
261          */
262         gmx_ana_poscalc_t *createCalculation(e_poscalc_t type, int flags);
263         /*! \brief
264          * Creates a new position calculation based on an enum value.
265          *
266          * \param[in]  post  One of the strings acceptable for
267          *      typeFromEnum().
268          * \param[in]  flags Flags for setting calculation options
269          *      (see \ref poscalc_flags "documentation of the flags").
270          * \throws     InternalError  if post is not recognized.
271          *
272          * This is a convenience wrapper for createCalculation().
273          * \p flags sets the default calculation options if not overridden by
274          * \p post; see typeFromEnum().
275          *
276          * May also throw std::bad_alloc in the future.
277          *
278          * \see createCalculation(), typeFromEnum()
279          */
280         gmx_ana_poscalc_t *createCalculationFromEnum(const char *post, int flags);
281
282         /*! \brief
283          * Computes the highest atom index required to evaluate this collection.
284          *
285          * Does not throw.
286          */
287         int getHighestRequiredAtomIndex() const;
288
289         /*! \brief
290          * Initializes evaluation for a position calculation collection.
291          *
292          * This function does some final initialization of the data structures
293          * in the collection to prepare them for evaluation.
294          * After this function has been called, it is no longer possible to add
295          * new calculations to the collection.
296          *
297          * Multiple calls to the function are ignored.
298          *
299          * Does not throw currently, but may throw std::bad_alloc in the
300          * future.
301          */
302         void initEvaluation();
303         /*! \brief
304          * Initializes a position calculation collection for a new frame.
305          *
306          * Clears the evaluation flag for all calculations.
307          * Should be called for each frame before calling
308          * gmx_ana_poscalc_update().
309          *
310          * This function calls initEvaluation() automatically if it has not
311          * been called earlier.
312          *
313          * Does not throw, but may throw if initEvaluation() is changed to
314          * throw.
315          */
316         void initFrame();
317
318     private:
319         class Impl;
320
321         PrivateImplPointer<Impl> impl_;
322
323         /*! \brief
324          * Needed to access the implementation class from the C code.
325          */
326         friend struct ::gmx_ana_poscalc_t;
327 };
328
329 } // namespace gmx
330
331 /** Sets the flags for position calculation. */
332 void
333 gmx_ana_poscalc_set_flags(gmx_ana_poscalc_t *pc, int flags);
334 /** Sets the maximum possible input index group for position calculation. */
335 void
336 gmx_ana_poscalc_set_maxindex(gmx_ana_poscalc_t *pc, gmx_ana_index_t *g);
337 /** Initializes positions for position calculation output. */
338 void
339 gmx_ana_poscalc_init_pos(gmx_ana_poscalc_t *pc, gmx_ana_pos_t *p);
340 /** Frees the memory allocated for position calculation. */
341 void
342 gmx_ana_poscalc_free(gmx_ana_poscalc_t *pc);
343 /** Returns true if the position calculation requires topology information. */
344 bool
345 gmx_ana_poscalc_requires_top(gmx_ana_poscalc_t *pc);
346
347 /** Updates a single COM/COG structure for a frame. */
348 void
349 gmx_ana_poscalc_update(gmx_ana_poscalc_t *pc,
350                        gmx_ana_pos_t *p, gmx_ana_index_t *g,
351                        t_trxframe *fr, t_pbc *pbc);
352
353 #endif