Make stepWorkload.useGpuXBufferOps flag consistent
[alexxy/gromacs.git] / src / gromacs / mdlib / vsite.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
11  *
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.
16  *
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.
21  *
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.
26  *
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.
34  *
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.
37  */
38 /*! \libinternal \file
39  * \brief Declares the VirtualSitesHandler class and vsite standalone functions
40  *
41  * \author Berk Hess <hess@kth.se>
42  * \ingroup module_mdlib
43  * \inlibraryapi
44  */
45
46 #ifndef GMX_MDLIB_VSITE_H
47 #define GMX_MDLIB_VSITE_H
48
49 #include <memory>
50
51 #include "gromacs/math/vectypes.h"
52 #include "gromacs/topology/idef.h"
53 #include "gromacs/utility/real.h"
54
55 struct gmx_domdec_t;
56 struct gmx_mtop_t;
57 struct t_commrec;
58 struct InteractionList;
59 struct t_nrnb;
60 struct gmx_wallcycle;
61 enum class PbcType : int;
62 enum class ParticleType : int;
63
64 namespace gmx
65 {
66 class RangePartitioning;
67 template<typename T>
68 class ArrayRef;
69
70 /*! \brief The start value of the vsite indices in the ftype enum
71  *
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.)
75  */
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;
79
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;
82
83 //! Whether we calculate vsite positions, velocities, or both
84 enum class VSiteOperation
85 {
86     Positions,              //!< Calculate only positions
87     Velocities,             //!< Calculate only velocities
88     PositionsAndVelocities, //!< Calculate both positions and velocities
89     Count                   //!< The number of entries
90 };
91
92 /*! \libinternal
93  * \brief Class that handles construction of vsites and spreading of vsite forces
94  */
95 class VirtualSitesHandler
96 {
97 public:
98     //! Constructor, used only be the makeVirtualSitesHandler() factory function
99     VirtualSitesHandler(const gmx_mtop_t&                 mtop,
100                         gmx_domdec_t*                     domdec,
101                         PbcType                           pbcType,
102                         ArrayRef<const RangePartitioning> updateGroupingPerMoleculeType);
103
104     ~VirtualSitesHandler();
105
106     //! Returns the number of virtual sites acting over multiple update groups
107     int numInterUpdategroupVirtualSites() const;
108
109     //! Set VSites and distribute VSite work over threads, should be called after each DD partitioning
110     void setVirtualSites(ArrayRef<const InteractionList> ilist,
111                          int                             numAtoms,
112                          int                             homenr,
113                          ArrayRef<const ParticleType>    ptype);
114
115     /*! \brief Create positions of vsite atoms based for the local system
116      *
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
121      */
122     void construct(ArrayRef<RVec> x, ArrayRef<RVec> v, const matrix box, VSiteOperation operation) const;
123
124     //! Tells how to handle virial contributions due to virtual sites
125     enum class VirialHandling : int
126     {
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
130     };
131
132     /*! \brief Spread the force operating on the vsite atoms on the surrounding atoms.
133      *
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.
141      */
142     void spreadForces(ArrayRef<const RVec> x,
143                       ArrayRef<RVec>       f,
144                       VirialHandling       virialHandling,
145                       ArrayRef<RVec>       fshift,
146                       matrix               virial,
147                       t_nrnb*              nrnb,
148                       const matrix         box,
149                       gmx_wallcycle*       wcycle);
150
151 private:
152     //! Implementation type.
153     class Impl;
154     //! Implementation object.
155     std::unique_ptr<Impl> impl_;
156 };
157
158 /*! \brief Create positions of vsite atoms based for the local system
159  *
160  * \param[in,out] x        The coordinates
161  * \param[in]     ip       Interaction parameters
162  * \param[in]     ilist    The interaction list
163  */
164 void constructVirtualSites(ArrayRef<RVec>                  x,
165                            ArrayRef<const t_iparams>       ip,
166                            ArrayRef<const InteractionList> ilist);
167
168 /*! \brief Create positions of vsite atoms for the whole system assuming all molecules are wholex
169  *
170  * \param[in]     mtop  The global topology
171  * \param[in,out] x     The global coordinates
172  */
173 void constructVirtualSitesGlobal(const gmx_mtop_t& mtop, ArrayRef<RVec> x);
174
175 //! Tells how to handle virial contributions due to virtual sites
176 enum class VirtualSiteVirialHandling : int
177 {
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
181 };
182
183 //! Return the number of non-linear virtual site constructions in the system
184 int countNonlinearVsites(const gmx_mtop_t& mtop);
185
186 /*! \brief Return the number of virtual sites that cross update groups
187  *
188  * \param[in] mtop                           The global topology
189  * \param[in] updateGroupingsPerMoleculeType  Update grouping per molecule type, pass empty when not using update groups
190  */
191 int countInterUpdategroupVsites(const gmx_mtop_t&                 mtop,
192                                 ArrayRef<const RangePartitioning> updateGroupingsPerMoleculeType);
193
194 /*! \brief Create the virtual site handler
195  *
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
202  */
203 std::unique_ptr<VirtualSitesHandler>
204 makeVirtualSitesHandler(const gmx_mtop_t&                 mtop,
205                         const t_commrec*                  cr,
206                         PbcType                           pbcType,
207                         ArrayRef<const RangePartitioning> updateGroupingPerMoleculeType);
208
209 } // namespace gmx
210
211 #endif