Use gmx_mtop_t in selections, part 2
[gromacs.git] / src / gromacs / selection / poscalc.h
blobe37628dea0fb28efeedde4accbb33cc37152853d
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2009,2010,2011,2012,2013,2014,2015,2016, 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.
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.
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.
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.
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.
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.
35 /*! \internal \file
36 * \brief
37 * API for structured and optimized calculation of positions.
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.
47 * The API is documented in more detail in gmx::PositionCalculationCollection.
49 * \author Teemu Murtola <teemu.murtola@gmail.com>
50 * \ingroup module_selection
52 #ifndef GMX_SELECTION_POSCALC_H
53 #define GMX_SELECTION_POSCALC_H
55 #include <cstdio>
57 #include "gromacs/utility/classhelpers.h"
59 /*! \name Flags for position calculation.
60 * \anchor poscalc_flags
62 /*@{*/
63 /*! \brief
64 * Use mass weighting.
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.
71 #define POS_MASS 1
72 /*! \brief
73 * Calculate positions for the same atoms in residues/molecules.
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.
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.
84 #define POS_COMPLMAX 2
85 /*! \brief
86 * Calculate positions for whole residues/molecules.
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.
92 * Has no effect unless the calculation type is \ref POS_RES or \ref POS_MOL.
94 #define POS_COMPLWHOLE 4
95 /*! \brief
96 * Enable handling of changing calculation groups.
98 * Can be used for static calculations as well, but implies a small
99 * performance penalty.
101 #define POS_DYNAMIC 16
102 /*! \brief
103 * Update \c gmx_ana_pos_t::m dynamically for an otherwise static
104 * calculation.
106 * Has effect only if \ref POS_DYNAMIC is not set.
108 #define POS_MASKONLY 32
109 /*! \brief
110 * Calculate velocities of the positions.
112 #define POS_VELOCITIES 64
113 /*! \brief
114 * Calculate forces on the positions.
116 #define POS_FORCES 128
117 /*@}*/
119 /** Specifies the type of positions to be calculated. */
120 typedef enum
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;
129 /** Data structure for position calculation. */
130 struct gmx_ana_poscalc_t;
132 struct gmx_ana_index_t;
133 struct gmx_ana_pos_t;
134 struct gmx_mtop_t;
135 struct t_pbc;
136 struct t_trxframe;
138 namespace gmx
141 /*! \internal
142 * \brief
143 * Collection of \c gmx_ana_poscalc_t structures for the same topology.
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.
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().
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.
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.
176 * \ingroup module_selection
178 class PositionCalculationCollection
180 public:
181 //! Describes what topology information is needed for position calculation.
182 enum class RequiredTopologyInfo
184 None, //!< No topology is needed.
185 Topology, //!< Topology is needed (residue/molecule info).
186 TopologyAndMasses //!< Masses are needed.
189 /*! \brief
190 * Array of strings acceptable for position calculation type enum.
192 * This array contains the acceptable values for typeFromEnum() and
193 * createCalculationFromEnum().
194 * The array contains a NULL pointer after the last item to indicate
195 * the end of the list.
197 static const char * const typeEnumValues[];
199 /*! \brief
200 * Converts a string to parameters for createCalculationFromEnum().
202 * \param[in] post String (typically an enum argument).
203 * Allowed values: 'atom', 'res_com', 'res_cog', 'mol_com', 'mol_cog',
204 * or one of the last four prepended by 'whole_', 'part_', or 'dyn_'.
205 * \param[out] type \c e_poscalc_t corresponding to \p post.
206 * \param[in,out] flags Flags corresponding to \p post.
207 * On input, the flags should contain the default flags.
208 * On exit, the flags \ref POS_MASS, \ref POS_COMPLMAX and
209 * \ref POS_COMPLWHOLE have been set according to \p post
210 * (the completion flags are left at the default values if no
211 * completion prefix is given).
212 * \throws InternalError if post is not recognized.
214 * \attention
215 * Checking is not complete, and other values than those listed above
216 * may be accepted for \p post, but the results are undefined.
218 * \see typeEnumValues
220 static void typeFromEnum(const char *post, e_poscalc_t *type, int *flags);
221 /*! \brief
222 * Returns what information is needed for position evaluation.
224 * \param[in] post Position type (see typeFromEnum()).
225 * \param[in] forces Whether forces are needed.
226 * \returns What topology information is required for initializing
227 * and/or evaluating the positions.
229 static RequiredTopologyInfo requiredTopologyInfoForType(const char *post, bool forces);
231 /*! \brief
232 * Creates a new position calculation collection object.
234 * \throws std::bad_alloc if out of memory.
236 PositionCalculationCollection();
237 /*! \brief
238 * Destroys a position calculation collection and its calculations.
240 * Any calculations in the collection are also freed, even if
241 * references to them are left.
243 ~PositionCalculationCollection();
245 /*! \brief
246 * Sets the topology used for the calculations.
248 * \param[in] top Topology data structure.
250 * This function should be called to set the topology before using
251 * gmx_ana_poscalc_set_maxindex() for any calculation that requires
252 * topology information.
254 * Does not throw.
256 void setTopology(const gmx_mtop_t *top);
257 /*! \brief
258 * Prints information about calculations.
260 * \param[in] fp File handle to receive the output.
262 * The output is very technical, making this function mainly useful for
263 * debugging purposes.
265 * Does not throw.
267 void printTree(FILE *fp) const;
269 /*! \brief
270 * Creates a new position calculation.
272 * \param[in] type Type of calculation.
273 * \param[in] flags Flags for setting calculation options
274 * (see \ref poscalc_flags "documentation of the flags").
276 * Does not throw currently, but may throw std::bad_alloc in the
277 * future.
279 gmx_ana_poscalc_t *createCalculation(e_poscalc_t type, int flags);
280 /*! \brief
281 * Creates a new position calculation based on an enum value.
283 * \param[in] post One of the strings acceptable for
284 * typeFromEnum().
285 * \param[in] flags Flags for setting calculation options
286 * (see \ref poscalc_flags "documentation of the flags").
287 * \throws InternalError if post is not recognized.
289 * This is a convenience wrapper for createCalculation().
290 * \p flags sets the default calculation options if not overridden by
291 * \p post; see typeFromEnum().
293 * May also throw std::bad_alloc in the future.
295 * \see createCalculation(), typeFromEnum()
297 gmx_ana_poscalc_t *createCalculationFromEnum(const char *post, int flags);
299 /*! \brief
300 * Computes the atoms required to evaluate this collection.
302 * \param[out] out Maximal group of atoms required to evaluate the
303 * positions.
305 * Does not throw.
307 void getRequiredAtoms(gmx_ana_index_t *out) const;
309 /*! \brief
310 * Initializes evaluation for a position calculation collection.
312 * This function does some final initialization of the data structures
313 * in the collection to prepare them for evaluation.
314 * After this function has been called, it is no longer possible to add
315 * new calculations to the collection.
317 * Multiple calls to the function are ignored.
319 * Does not throw currently, but may throw std::bad_alloc in the
320 * future.
322 void initEvaluation();
323 /*! \brief
324 * Initializes a position calculation collection for a new frame.
326 * \param[in] fr Frame to initialize evaluation for.
328 * Should be called for each frame before calling
329 * gmx_ana_poscalc_update().
331 * This function calls initEvaluation() automatically if it has not
332 * been called earlier.
334 * Does not currently throw, but may throw std::bad_alloc in the
335 * future.
337 void initFrame(const t_trxframe *fr);
339 private:
340 class Impl;
342 PrivateImplPointer<Impl> impl_;
344 /*! \brief
345 * Needed to access the implementation class from the C code.
347 friend struct ::gmx_ana_poscalc_t;
350 } // namespace gmx
352 /** Sets the flags for position calculation. */
353 void
354 gmx_ana_poscalc_set_flags(gmx_ana_poscalc_t *pc, int flags);
355 /** Sets the maximum possible input index group for position calculation. */
356 void
357 gmx_ana_poscalc_set_maxindex(gmx_ana_poscalc_t *pc, gmx_ana_index_t *g);
358 /** Initializes positions for position calculation output. */
359 void
360 gmx_ana_poscalc_init_pos(gmx_ana_poscalc_t *pc, gmx_ana_pos_t *p);
361 /** Frees the memory allocated for position calculation. */
362 void
363 gmx_ana_poscalc_free(gmx_ana_poscalc_t *pc);
364 /*! \brief
365 * Returns true if the position calculation requires topology information.
367 * \param[in] pc Position calculation data to query.
368 * \returns Which topology information \p pc requires for initialization
369 * and/or evaluation.
371 gmx::PositionCalculationCollection::RequiredTopologyInfo
372 gmx_ana_poscalc_required_topology_info(gmx_ana_poscalc_t *pc);
374 /** Updates a single COM/COG structure for a frame. */
375 void
376 gmx_ana_poscalc_update(gmx_ana_poscalc_t *pc,
377 gmx_ana_pos_t *p, gmx_ana_index_t *g,
378 t_trxframe *fr, const t_pbc *pbc);
380 #endif