Split simulationWork.useGpuBufferOps into separate x and f flags
[alexxy/gromacs.git] / src / gromacs / mdlib / settle_gpu.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2019,2020,2021, 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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 /*! \internal \file
36  *
37  * \brief Declares class for GPU implementation of SETTLE
38  *
39  * \author Artem Zhmurov <zhmurov@gmail.com>
40  *
41  * \ingroup module_mdlib
42  */
43 #ifndef GMX_MDLIB_SETTLE_GPU_H
44 #define GMX_MDLIB_SETTLE_GPU_H
45
46 #include "gmxpre.h"
47
48 #include "gromacs/gpu_utils/devicebuffer_datatype.h"
49 #include "gromacs/gpu_utils/device_context.h"
50 #include "gromacs/gpu_utils/device_stream.h"
51 #include "gromacs/gpu_utils/gputraits.h"
52 #include "gromacs/math/functions.h"
53 #include "gromacs/math/invertmatrix.h"
54 #include "gromacs/math/vec.h"
55 #include "gromacs/mdlib/settle.h"
56 #include "gromacs/pbcutil/pbc.h"
57 #include "gromacs/pbcutil/pbc_aiuc.h"
58 #include "gromacs/topology/topology.h"
59
60 class InteractionDefinitions;
61
62 namespace gmx
63 {
64
65 //! \internal \brief Indices of atoms in a water molecule
66 struct WaterMolecule
67 {
68     //! Oxygen atom
69     int ow1;
70     //! First hydrogen atom
71     int hw2;
72     //! Second hydrogen atom
73     int hw3;
74 };
75
76
77 /*! \internal \brief Class with interfaces and data for GPU version of SETTLE. */
78 class SettleGpu
79 {
80
81 public:
82     /*! \brief Create SETTLE object
83      *
84      *  Extracts masses for oxygen and hydrogen as well as the O-H and H-H target distances
85      *  from the topology data (mtop), check their values for consistency and calls the
86      *  following constructor.
87      *
88      * \param[in] mtop           Topology of the system to gen the masses for O and H atoms and
89      *                           target O-H and H-H distances. These values are also checked for
90      *                           consistency.
91      * \param[in] deviceContext  Device context (dummy in CUDA).
92      * \param[in] deviceStream   Device stream to use.
93      */
94     SettleGpu(const gmx_mtop_t& mtop, const DeviceContext& deviceContext, const DeviceStream& deviceStream);
95
96     ~SettleGpu();
97
98     /*! \brief Apply SETTLE.
99      *
100      * Applies SETTLE to coordinates and velocities, stored on GPU. Data at pointers d_xp and
101      * d_v change in the GPU memory. The results are not automatically copied back to the CPU
102      * memory. Method uses this class data structures which should be updated when needed using
103      * update method.
104      *
105      * \param[in]     d_x               Coordinates before timestep (in GPU memory)
106      * \param[in,out] d_xp              Coordinates after timestep (in GPU memory). The
107      *                                  resulting constrained coordinates will be saved here.
108      * \param[in]     updateVelocities  If the velocities should be updated.
109      * \param[in,out] d_v               Velocities to update (in GPU memory, can be nullptr
110      *                                  if not updated)
111      * \param[in]     invdt             Reciprocal timestep (to scale Lagrange
112      *                                  multipliers when velocities are updated)
113      * \param[in]     computeVirial     If virial should be updated.
114      * \param[in,out] virialScaled      Scaled virial tensor to be updated.
115      * \param[in]     pbcAiuc           PBC data.
116      */
117     void apply(const DeviceBuffer<Float3>& d_x,
118                DeviceBuffer<Float3>        d_xp,
119                bool                        updateVelocities,
120                DeviceBuffer<Float3>        d_v,
121                real                        invdt,
122                bool                        computeVirial,
123                tensor                      virialScaled,
124                const PbcAiuc&              pbcAiuc);
125
126     /*! \brief
127      * Update data-structures (e.g. after NB search step).
128      *
129      * Updates the constraints data and copies it to the GPU. Should be
130      * called if the particles were sorted, redistributed between domains, etc.
131      * Does not recycle the data preparation routines from the CPU version.
132      * All three atoms from single water molecule should be handled by the same GPU.
133      *
134      * SETTLEs atom ID's is taken from idef.il[F_SETTLE].iatoms.
135      *
136      * \param[in] idef    System topology
137      */
138     void set(const InteractionDefinitions& idef);
139
140 private:
141     //! GPU context object
142     const DeviceContext& deviceContext_;
143     //! GPU stream
144     const DeviceStream& deviceStream_;
145
146     //! Scaled virial tensor (9 reals, GPU)
147     std::vector<float> h_virialScaled_;
148     //! Scaled virial tensor (9 reals, GPU)
149     DeviceBuffer<float> d_virialScaled_;
150
151     //! Number of settles
152     int numSettles_ = 0;
153
154     //! Indexes of atoms (.i for oxygen, .j and.k for hydrogens, CPU)
155     std::vector<WaterMolecule> h_atomIds_;
156     //! Indexes of atoms (.i for oxygen, .j and.k for hydrogens, GPU)
157     DeviceBuffer<WaterMolecule> d_atomIds_;
158     //! Current size of the array of atom IDs
159     int numAtomIds_ = -1;
160     //! Allocated size for the array of atom IDs
161     int numAtomIdsAlloc_ = -1;
162
163     //! Settle parameters
164     SettleParameters settleParameters_;
165 };
166
167 } // namespace gmx
168
169 #endif // GMX_MDLIB_SETTLE_GPU_H