2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2018,2019, 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 /*! \libinternal \file
37 * \brief Declares the UpdateGroupsCog class for managing centers of mass of update groups
39 * \author Berk Hess <hess@kth.se>
40 * \ingroup module_mdlib
43 #ifndef GMX_MDLIB_UPDATEGROUPSCOG
44 #define GMX_MDLIB_UPDATEGROUPSCOG
48 #include "gromacs/domdec/hashedmap.h"
49 #include "gromacs/math/vectypes.h"
50 #include "gromacs/utility/arrayref.h"
51 #include "gromacs/utility/gmxassert.h"
56 class RangePartitioning;
59 * \brief Class for managing and computing centers of geometry of update groups
64 /*! \brief Constructor
66 * The temperature is used for computing the maximum update group
69 * \note \p numHomeAtoms only affects the performance up till the first
72 * \param[in] mtop The global topology
73 * \param[in] updateGroupsPerMoleculetype List of update groups for each molecule type in \p mtop
74 * \param[in] temperature The maximum reference temperature, pass -1 when unknown or not applicable
75 * \param[in] numHomeAtoms Estimate of the number of home atoms per DD cell
77 UpdateGroupsCog(const gmx_mtop_t& mtop,
78 gmx::ArrayRef<const gmx::RangePartitioning> updateGroupsPerMoleculetype,
82 /*! \brief Compute centers of geometry for supplied coordinates
84 * Coordinates are processed starting after the last index
85 * processed in the previous call to \p addCogs(), unless \p clear()
86 * was called last, in which case processing starts at 0.
88 * \param[in] globalAtomIndices List of global atom indices for the atoms belonging to \p coordinates
89 * \param[in] coordinates List of coordinates to be processed, processing might not start at 0 (see above)
91 void addCogs(gmx::ArrayRef<const int> globalAtomIndices, gmx::ArrayRef<const gmx::RVec> coordinates);
93 /*! \brief Returns the number of centers of geometry currently stored */
94 int numCogs() const { return cogs_.size(); }
96 /*! \brief Returns a reference to a center of geometry
98 * \param[in] cogIndex The COG requested, should be; 0 <= \p cogIndex < cogs_.size()
100 RVec& cog(int cogIndex)
102 GMX_ASSERT(cogIndex >= 0 && static_cast<size_t>(cogIndex) < cogs_.size(),
103 "cogIndex should be in the range set in this object");
105 return cogs_[cogIndex];
108 /*! \brief Returns the COG index given an atom index
110 * \param[in] atomIndex The local atom index that maps to the COG
112 int cogIndex(int atomIndex) const
114 GMX_ASSERT(atomIndex >= 0 && static_cast<size_t>(atomIndex) < cogIndices_.size(),
115 "atomIndex should be in the range set in this object");
117 return cogIndices_[atomIndex];
120 /*! \brief Return a reference to a center of geometry given an atom index
122 * \param[in] atomIndex The atom index that maps to the COG requested
124 const RVec& cogForAtom(int atomIndex) const { return cogs_[cogIndex(atomIndex)]; }
126 /*! \brief Clear the lists of stored center of geometry coordinates */
129 /*! \brief Returns the maximum radius over all update groups */
130 real maxUpdateGroupRadius() const { return maxUpdateGroupRadius_; }
134 * \brief Helper struct for mapping atom indices to COG indices, used per molecule block in mtop
138 //! \brief Starting index of groups for this molblock
140 //! \brief The number of atoms in a molecule
141 int numGroupsPerMolecule_;
142 //! \brief Map atom indices to group indices for one molecule
143 std::vector<int> groupIndex_;
146 //! \brief Maps atom indices to COG indices
147 std::vector<int> cogIndices_;
148 //! \brief List of COGs
149 std::vector<RVec> cogs_;
150 //! \brief List of the number of atoms for each COG
151 std::vector<int> numAtomsPerCog_;
152 //! \brief Maps global COG index to local COG index
153 HashedMap<int> globalToLocalMap_;
154 //! \brief Helper data for mapping atom indices to COG indices
155 std::vector<IndexToGroup> indicesPerMoleculeblock_;
156 //! \brief The maximum radius over all update groups
157 real maxUpdateGroupRadius_;
158 //! \brief Reference to mtop this object was constructed with
159 const gmx_mtop_t& mtop_;
164 #endif // GMX_MDLIB_UPDATEGROUPSCOG