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/arrayref.h"
54 #include "gromacs/utility/basedefinitions.h"
55 #include "gromacs/utility/real.h"
60 struct InteractionList;
64 enum class PbcType : int;
68 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, gmx_domdec_t* domdec, PbcType pbcType);
101 ~VirtualSitesHandler();
103 //! Returns the number of virtual sites acting over multiple update groups
104 int numInterUpdategroupVirtualSites() const;
106 //! Set VSites and distribute VSite work over threads, should be called after each DD partitioning
107 void setVirtualSites(ArrayRef<const InteractionList> ilist, const t_mdatoms& mdatoms);
109 /*! \brief Create positions of vsite atoms based for the local system
111 * \param[in,out] x The coordinates
112 * \param[in,out] v The velocities, needed if operation requires it
113 * \param[in] box The box
114 * \param[in] operation Whether we calculate positions, velocities, or both
116 void construct(ArrayRef<RVec> x, ArrayRef<RVec> v, const matrix box, VSiteOperation operation) const;
118 //! Tells how to handle virial contributions due to virtual sites
119 enum class VirialHandling : int
121 None, //!< Do not compute virial contributions
122 Pbc, //!< Add contributions working over PBC to shift forces
123 NonLinear //!< Compute contributions due to non-linear virtual sites
126 /*! \brief Spread the force operating on the vsite atoms on the surrounding atoms.
128 * vsite should point to a valid object.
129 * The virialHandling parameter determines how virial contributions are handled.
130 * If this is set to Linear, shift forces are accumulated into fshift.
131 * If this is set to NonLinear, non-linear contributions are added to virial.
132 * This non-linear correction is required when the virial is not calculated
133 * afterwards from the particle position and forces, but in a different way,
134 * as for instance for the PME mesh contribution.
136 void spreadForces(ArrayRef<const RVec> x,
138 VirialHandling virialHandling,
139 ArrayRef<RVec> fshift,
143 gmx_wallcycle* wcycle);
146 //! Implementation type.
148 //! Implementation object.
149 std::unique_ptr<Impl> impl_;
152 /*! \brief Create positions of vsite atoms based for the local system
154 * \param[in,out] x The coordinates
155 * \param[in] ip Interaction parameters
156 * \param[in] ilist The interaction list
158 void constructVirtualSites(ArrayRef<RVec> x,
159 ArrayRef<const t_iparams> ip,
160 ArrayRef<const InteractionList> ilist);
162 /*! \brief Create positions of vsite atoms for the whole system assuming all molecules are wholex
164 * \param[in] mtop The global topology
165 * \param[in,out] x The global coordinates
167 void constructVirtualSitesGlobal(const gmx_mtop_t& mtop, ArrayRef<RVec> x);
169 //! Tells how to handle virial contributions due to virtual sites
170 enum class VirtualSiteVirialHandling : int
172 None, //!< Do not compute virial contributions
173 Pbc, //!< Add contributions working over PBC to shift forces
174 NonLinear //!< Compute contributions due to non-linear virtual sites
177 //! Return the number of non-linear virtual site constructions in the system
178 int countNonlinearVsites(const gmx_mtop_t& mtop);
180 /*! \brief Return the number of virtual sites that cross update groups
182 * \param[in] mtop The global topology
183 * \param[in] updateGroupingPerMoleculetype Update grouping per molecule type, pass empty when not using update groups
185 int countInterUpdategroupVsites(const gmx_mtop_t& mtop,
186 ArrayRef<const RangePartitioning> updateGroupingPerMoleculetype);
188 /*! \brief Create the virtual site handler
190 * \param[in] mtop The global topology
191 * \param[in] cr The communication record
192 * \param[in] pbcType The type of PBC
193 * \returns A valid vsite handler object or nullptr when there are no virtual sites
195 std::unique_ptr<VirtualSitesHandler> makeVirtualSitesHandler(const gmx_mtop_t& mtop,