Apply clang-format to source tree
[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-2016, The GROMACS development team.
5  * Copyright (c) 2019, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36 /*! \internal \file
37  * \brief
38  * API for structured and optimized calculation of positions.
39  *
40  * This header declares an API for calculating positions in an automated way,
41  * for internal use by the selection engine.  This is useful in particular with
42  * dynamic selections, because the same COM/COG positions may be needed in
43  * several contexts.  The API makes it possible to optimize the evaluation such
44  * that any heavy calculation is only done once, and the results just copied if
45  * needed more than once.  The functions also provide a convenient interface
46  * for keeping the whole \c gmx_ana_pos_t structure up-to-date.
47  *
48  * The API is documented in more detail in gmx::PositionCalculationCollection.
49  *
50  * \author Teemu Murtola <teemu.murtola@gmail.com>
51  * \ingroup module_selection
52  */
53 #ifndef GMX_SELECTION_POSCALC_H
54 #define GMX_SELECTION_POSCALC_H
55
56 #include <cstdio>
57
58 #include "gromacs/utility/classhelpers.h"
59
60 /*! \name Flags for position calculation.
61  * \anchor poscalc_flags
62  */
63 /*@{*/
64 /*! \brief
65  * Use mass weighting.
66  *
67  * If this flag is set, the positions will be calculated using mass weighting,
68  * i.e., one gets center-of-mass positions.
69  * Without the flag, center-of-geometry positions are calculated.
70  * Does not have any effect if the calculation type is \ref POS_ATOM.
71  */
72 #define POS_MASS 1
73 /*! \brief
74  * Calculate positions for the same atoms in residues/molecules.
75  *
76  * If this flag is set, the positions are always calculated using the same
77  * atoms for each residue/molecule, even if the evaluation group contains only
78  * some of the atoms for some frames.
79  * The group passed to gmx_ana_poscalc_set_maxindex() is used to determine
80  * the atoms to use for the calculation.
81  *
82  * Has no effect unless \ref POS_DYNAMIC is set or if the calculation type
83  * is not \ref POS_RES of \ref POS_MOL.
84  */
85 #define POS_COMPLMAX 2
86 /*! \brief
87  * Calculate positions for whole residues/molecules.
88  *
89  * If this flag is set, the positions will be calculated for whole
90  * residues/molecules, even if the group contains only some of the atoms in
91  * the residue/molecule.
92  *
93  * Has no effect unless the calculation type is \ref POS_RES or \ref POS_MOL.
94  */
95 #define POS_COMPLWHOLE 4
96 /*! \brief
97  * Enable handling of changing calculation groups.
98  *
99  * Can be used for static calculations as well, but implies a small
100  * performance penalty.
101  */
102 #define POS_DYNAMIC 16
103 /*! \brief
104  * Update \c gmx_ana_pos_t::m dynamically for an otherwise static
105  * calculation.
106  *
107  * Has effect only if \ref POS_DYNAMIC is not set.
108  */
109 #define POS_MASKONLY 32
110 /*! \brief
111  * Calculate velocities of the positions.
112  */
113 #define POS_VELOCITIES 64
114 /*! \brief
115  * Calculate forces on the positions.
116  */
117 #define POS_FORCES 128
118 /*@}*/
119
120 /** Specifies the type of positions to be calculated. */
121 typedef enum
122 {
123     POS_ATOM,   /**< Copy atomic coordinates. */
124     POS_RES,    /**< Calculate center for each residue. */
125     POS_MOL,    /**< Calculate center for each molecule. */
126     POS_ALL,    /**< Calculate center for the whole group. */
127     POS_ALL_PBC /**< Calculate center for the whole group with PBC. */
128 } e_poscalc_t;
129
130 /** Data structure for position calculation. */
131 struct gmx_ana_poscalc_t;
132
133 struct gmx_ana_index_t;
134 struct gmx_ana_pos_t;
135 struct gmx_mtop_t;
136 struct t_pbc;
137 struct t_trxframe;
138
139 namespace gmx
140 {
141
142 /*! \internal
143  * \brief
144  * Collection of \c gmx_ana_poscalc_t structures for the same topology.
145  *
146  * Calculations within one collection share the same topology, and they are
147  * optimized.  Calculations in different collections do not interact.
148  * The topology for a collection can be set with setTopology().
149  * This needs to be done before calling gmx_ana_poscalc_set_maxindex() for
150  * any calculation in the collection, unless that calculation does not
151  * require topology information.
152  *
153  * A new calculation is created with createCalculation().
154  * If flags need to be adjusted later, gmx_ana_poscalc_set_flags() can be
155  * used.
156  * After the flags are final, the largest possible index group for which the
157  * positions are needed has to be set with gmx_ana_poscalc_set_maxindex().
158  * setTopology() should have been called before this function is called.
159  * After the above calls, gmx_ana_poscalc_init_pos() can be used to initialize
160  * output to a \c gmx_ana_pos_t structure.  Several different structures can be
161  * initialized for the same calculation; the only requirement is that the
162  * structure passed later to gmx_ana_poscalc_update() has been initialized
163  * properly.
164  * The memory allocated for a calculation can be freed with
165  * gmx_ana_poscalc_free().
166  *
167  * The position evaluation is simple: initFrame() should be
168  * called once for each frame, and gmx_ana_poscalc_update() can then be called
169  * for each calculation that is needed for that frame.
170  *
171  * It is also possible to initialize the calculations based on a type provided
172  * as a string.
173  * The possible strings are listed in \ref typeEnumValues, and the string can
174  * be converted to the parameters for createCalculation() using typeFromEnum().
175  * createCalculationFromEnum() is also provided for convenience.
176  *
177  * \ingroup module_selection
178  */
179 class PositionCalculationCollection
180 {
181 public:
182     //! Describes what topology information is needed for position calculation.
183     enum class RequiredTopologyInfo
184     {
185         None,             //!< No topology is needed.
186         Topology,         //!< Topology is needed (residue/molecule info).
187         TopologyAndMasses //!< Masses are needed.
188     };
189
190     /*! \brief
191      * Array of strings acceptable for position calculation type enum.
192      *
193      * This array contains the acceptable values for typeFromEnum() and
194      * createCalculationFromEnum().
195      * The array contains a NULL pointer after the last item to indicate
196      * the end of the list.
197      */
198     static const char* const typeEnumValues[];
199
200     /*! \brief
201      * Converts a string to parameters for createCalculationFromEnum().
202      *
203      * \param[in]     post  String (typically an enum argument).
204      *     Allowed values: 'atom', 'res_com', 'res_cog', 'mol_com', 'mol_cog',
205      *     or one of the last four prepended by 'whole_', 'part_', or 'dyn_'.
206      * \param[out]    type  \c e_poscalc_t corresponding to \p post.
207      * \param[in,out] flags Flags corresponding to \p post.
208      *     On input, the flags should contain the default flags.
209      *     On exit, the flags \ref POS_MASS, \ref POS_COMPLMAX and
210      *     \ref POS_COMPLWHOLE have been set according to \p post
211      *     (the completion flags are left at the default values if no
212      *     completion prefix is given).
213      * \throws  InternalError  if post is not recognized.
214      *
215      * \attention
216      * Checking is not complete, and other values than those listed above
217      * may be accepted for \p post, but the results are undefined.
218      *
219      * \see typeEnumValues
220      */
221     static void typeFromEnum(const char* post, e_poscalc_t* type, int* flags);
222     /*! \brief
223      * Returns what information is needed for position evaluation.
224      *
225      * \param[in] post   Position type (see typeFromEnum()).
226      * \param[in] forces Whether forces are needed.
227      * \returns   What topology information is required for initializing
228      *     and/or evaluating the positions.
229      */
230     static RequiredTopologyInfo requiredTopologyInfoForType(const char* post, bool forces);
231
232     /*! \brief
233      * Creates a new position calculation collection object.
234      *
235      * \throws  std::bad_alloc if out of memory.
236      */
237     PositionCalculationCollection();
238     /*! \brief
239      * Destroys a position calculation collection and its calculations.
240      *
241      * Any calculations in the collection are also freed, even if
242      * references to them are left.
243      */
244     ~PositionCalculationCollection();
245
246     /*! \brief
247      * Sets the topology used for the calculations.
248      *
249      * \param[in]     top   Topology data structure.
250      *
251      * This function should be called to set the topology before using
252      * gmx_ana_poscalc_set_maxindex() for any calculation that requires
253      * topology information.
254      *
255      * Does not throw.
256      */
257     void setTopology(const gmx_mtop_t* top);
258     /*! \brief
259      * Prints information about calculations.
260      *
261      * \param[in] fp    File handle to receive the output.
262      *
263      * The output is very technical, making this function mainly useful for
264      * debugging purposes.
265      *
266      * Does not throw.
267      */
268     void printTree(FILE* fp) const;
269
270     /*! \brief
271      * Creates a new position calculation.
272      *
273      * \param[in]  type  Type of calculation.
274      * \param[in]  flags Flags for setting calculation options
275      *   (see \ref poscalc_flags "documentation of the flags").
276      *
277      * Does not throw currently, but may throw std::bad_alloc in the
278      * future.
279      */
280     gmx_ana_poscalc_t* createCalculation(e_poscalc_t type, int flags);
281     /*! \brief
282      * Creates a new position calculation based on an enum value.
283      *
284      * \param[in]  post  One of the strings acceptable for
285      *      typeFromEnum().
286      * \param[in]  flags Flags for setting calculation options
287      *      (see \ref poscalc_flags "documentation of the flags").
288      * \throws     InternalError  if post is not recognized.
289      *
290      * This is a convenience wrapper for createCalculation().
291      * \p flags sets the default calculation options if not overridden by
292      * \p post; see typeFromEnum().
293      *
294      * May also throw std::bad_alloc in the future.
295      *
296      * \see createCalculation(), typeFromEnum()
297      */
298     gmx_ana_poscalc_t* createCalculationFromEnum(const char* post, int flags);
299
300     /*! \brief
301      * Computes the atoms required to evaluate this collection.
302      *
303      * \param[out] out  Maximal group of atoms required to evaluate the
304      *     positions.
305      *
306      * Does not throw.
307      */
308     void getRequiredAtoms(gmx_ana_index_t* out) const;
309
310     /*! \brief
311      * Initializes evaluation for a position calculation collection.
312      *
313      * This function does some final initialization of the data structures
314      * in the collection to prepare them for evaluation.
315      * After this function has been called, it is no longer possible to add
316      * new calculations to the collection.
317      *
318      * Multiple calls to the function are ignored.
319      *
320      * Does not throw currently, but may throw std::bad_alloc in the
321      * future.
322      */
323     void initEvaluation();
324     /*! \brief
325      * Initializes a position calculation collection for a new frame.
326      *
327      * \param[in] fr  Frame to initialize evaluation for.
328      *
329      * Should be called for each frame before calling
330      * gmx_ana_poscalc_update().
331      *
332      * This function calls initEvaluation() automatically if it has not
333      * been called earlier.
334      *
335      * Does not currently throw, but may throw std::bad_alloc in the
336      * future.
337      */
338     void initFrame(const t_trxframe* fr);
339
340 private:
341     class Impl;
342
343     PrivateImplPointer<Impl> impl_;
344
345     /*! \brief
346      * Needed to access the implementation class from the C code.
347      */
348     friend struct ::gmx_ana_poscalc_t;
349 };
350
351 } // namespace gmx
352
353 /** Sets the flags for position calculation. */
354 void 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 gmx_ana_poscalc_set_maxindex(gmx_ana_poscalc_t* pc, gmx_ana_index_t* g);
357 /** Initializes positions for position calculation output. */
358 void gmx_ana_poscalc_init_pos(gmx_ana_poscalc_t* pc, gmx_ana_pos_t* p);
359 /** Frees the memory allocated for position calculation. */
360 void gmx_ana_poscalc_free(gmx_ana_poscalc_t* pc);
361 /*! \brief
362  * Returns true if the position calculation requires topology information.
363  *
364  * \param[in] pc  Position calculation data to query.
365  * \returns   Which topology information \p pc requires for initialization
366  *     and/or evaluation.
367  */
368 gmx::PositionCalculationCollection::RequiredTopologyInfo
369 gmx_ana_poscalc_required_topology_info(gmx_ana_poscalc_t* pc);
370
371 /** Updates a single COM/COG structure for a frame. */
372 void gmx_ana_poscalc_update(gmx_ana_poscalc_t* pc,
373                             gmx_ana_pos_t*     p,
374                             gmx_ana_index_t*   g,
375                             t_trxframe*        fr,
376                             const t_pbc*       pbc);
377
378 #endif