2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016,2017 The GROMACS development team.
7 * Copyright (c) 2018,2019,2020,2021, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
38 /*! \libinternal \file
39 * \brief Declares the VirtualSitesHandler class and vsite standalone functions
41 * \author Berk Hess <hess@kth.se>
42 * \ingroup module_mdlib
46 #ifndef GMX_MDLIB_VSITE_H
47 #define GMX_MDLIB_VSITE_H
51 #include "gromacs/math/vectypes.h"
52 #include "gromacs/topology/idef.h"
53 #include "gromacs/utility/real.h"
58 struct InteractionList;
61 enum class PbcType : int;
62 enum class ParticleType : int;
66 class RangePartitioning;
70 /*! \brief The start value of the vsite indices in the ftype enum
72 * The validity of the start and end values is checked in makeVirtualSitesHandler().
73 * This is used to avoid loops over all ftypes just to get the vsite entries.
74 * (We should replace the fixed ilist array by only the used entries.)
76 static constexpr int c_ftypeVsiteStart = F_VSITE1;
77 //! The start and end value of the vsite indices in the ftype enum
78 static constexpr int c_ftypeVsiteEnd = F_VSITEN + 1;
80 //! Type for storing PBC atom information for all vsite types in the system
81 typedef std::array<std::vector<int>, c_ftypeVsiteEnd - c_ftypeVsiteStart> VsitePbc;
83 //! Whether we calculate vsite positions, velocities, or both
84 enum class VSiteOperation
86 Positions, //!< Calculate only positions
87 Velocities, //!< Calculate only velocities
88 PositionsAndVelocities, //!< Calculate both positions and velocities
89 Count //!< The number of entries
93 * \brief Class that handles construction of vsites and spreading of vsite forces
95 class VirtualSitesHandler
98 //! Constructor, used only be the makeVirtualSitesHandler() factory function
99 VirtualSitesHandler(const gmx_mtop_t& mtop,
100 gmx_domdec_t* domdec,
102 ArrayRef<const RangePartitioning> updateGroupingPerMoleculeType);
104 ~VirtualSitesHandler();
106 //! Returns the number of virtual sites acting over multiple update groups
107 int numInterUpdategroupVirtualSites() const;
109 //! Set VSites and distribute VSite work over threads, should be called after each DD partitioning
110 void setVirtualSites(ArrayRef<const InteractionList> ilist,
113 ArrayRef<const ParticleType> ptype);
115 /*! \brief Create positions of vsite atoms based for the local system
117 * \param[in,out] x The coordinates
118 * \param[in,out] v The velocities, needed if operation requires it
119 * \param[in] box The box
120 * \param[in] operation Whether we calculate positions, velocities, or both
122 void construct(ArrayRef<RVec> x, ArrayRef<RVec> v, const matrix box, VSiteOperation operation) const;
124 //! Tells how to handle virial contributions due to virtual sites
125 enum class VirialHandling : int
127 None, //!< Do not compute virial contributions
128 Pbc, //!< Add contributions working over PBC to shift forces
129 NonLinear //!< Compute contributions due to non-linear virtual sites
132 /*! \brief Spread the force operating on the vsite atoms on the surrounding atoms.
134 * vsite should point to a valid object.
135 * The virialHandling parameter determines how virial contributions are handled.
136 * If this is set to Linear, shift forces are accumulated into fshift.
137 * If this is set to NonLinear, non-linear contributions are added to virial.
138 * This non-linear correction is required when the virial is not calculated
139 * afterwards from the particle position and forces, but in a different way,
140 * as for instance for the PME mesh contribution.
142 void spreadForces(ArrayRef<const RVec> x,
144 VirialHandling virialHandling,
145 ArrayRef<RVec> fshift,
149 gmx_wallcycle* wcycle);
152 //! Implementation type.
154 //! Implementation object.
155 std::unique_ptr<Impl> impl_;
158 /*! \brief Create positions of vsite atoms based for the local system
160 * \param[in,out] x The coordinates
161 * \param[in] ip Interaction parameters
162 * \param[in] ilist The interaction list
164 void constructVirtualSites(ArrayRef<RVec> x,
165 ArrayRef<const t_iparams> ip,
166 ArrayRef<const InteractionList> ilist);
168 /*! \brief Create positions of vsite atoms for the whole system assuming all molecules are wholex
170 * \param[in] mtop The global topology
171 * \param[in,out] x The global coordinates
173 void constructVirtualSitesGlobal(const gmx_mtop_t& mtop, ArrayRef<RVec> x);
175 //! Tells how to handle virial contributions due to virtual sites
176 enum class VirtualSiteVirialHandling : int
178 None, //!< Do not compute virial contributions
179 Pbc, //!< Add contributions working over PBC to shift forces
180 NonLinear //!< Compute contributions due to non-linear virtual sites
183 //! Return the number of non-linear virtual site constructions in the system
184 int countNonlinearVsites(const gmx_mtop_t& mtop);
186 /*! \brief Return the number of virtual sites that cross update groups
188 * \param[in] mtop The global topology
189 * \param[in] updateGroupingsPerMoleculeType Update grouping per molecule type, pass empty when not using update groups
191 int countInterUpdategroupVsites(const gmx_mtop_t& mtop,
192 ArrayRef<const RangePartitioning> updateGroupingsPerMoleculeType);
194 /*! \brief Create the virtual site handler
196 * \param[in] mtop The global topology
197 * \param[in] cr The communication record
198 * \param[in] pbcType The type of PBC
199 * \param[in] updateGroupingPerMoleculeType Update grouping per molecule type, pass
200 * empty when not using update groups
201 * \returns A valid vsite handler object or nullptr when there are no virtual sites
203 std::unique_ptr<VirtualSitesHandler>
204 makeVirtualSitesHandler(const gmx_mtop_t& mtop,
207 ArrayRef<const RangePartitioning> updateGroupingPerMoleculeType);