f09fbd344f4a37b96c5e21fc7fcebade243fd0c4
[alexxy/gromacs.git] / src / gromacs / mdlib / settle_gpu.cuh
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_CUH
44 #define GMX_MDLIB_SETTLE_GPU_CUH
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 Class with interfaces and data for GPU version of SETTLE. */
66 class SettleGpu
67 {
68
69 public:
70     /*! \brief Create SETTLE object
71      *
72      *  Extracts masses for oxygen and hydrogen as well as the O-H and H-H target distances
73      *  from the topology data (mtop), check their values for consistency and calls the
74      *  following constructor.
75      *
76      * \param[in] mtop           Topology of the system to gen the masses for O and H atoms and
77      *                           target O-H and H-H distances. These values are also checked for
78      *                           consistency.
79      * \param[in] deviceContext  Device context (dummy in CUDA).
80      * \param[in] deviceStream   Device stream to use.
81      */
82     SettleGpu(const gmx_mtop_t& mtop, const DeviceContext& deviceContext, const DeviceStream& deviceStream);
83
84     ~SettleGpu();
85
86     /*! \brief Apply SETTLE.
87      *
88      * Applies SETTLE to coordinates and velocities, stored on GPU. Data at pointers d_xp and
89      * d_v change in the GPU memory. The results are not automatically copied back to the CPU
90      * memory. Method uses this class data structures which should be updated when needed using
91      * update method.
92      *
93      * \param[in]     d_x               Coordinates before timestep (in GPU memory)
94      * \param[in,out] d_xp              Coordinates after timestep (in GPU memory). The
95      *                                  resulting constrained coordinates will be saved here.
96      * \param[in]     updateVelocities  If the velocities should be updated.
97      * \param[in,out] d_v               Velocities to update (in GPU memory, can be nullptr
98      *                                  if not updated)
99      * \param[in]     invdt             Reciprocal timestep (to scale Lagrange
100      *                                  multipliers when velocities are updated)
101      * \param[in]     computeVirial     If virial should be updated.
102      * \param[in,out] virialScaled      Scaled virial tensor to be updated.
103      * \param[in]     pbcAiuc           PBC data.
104      */
105     void apply(const DeviceBuffer<Float3> d_x,
106                DeviceBuffer<Float3>       d_xp,
107                const bool                 updateVelocities,
108                DeviceBuffer<Float3>       d_v,
109                const real                 invdt,
110                const bool                 computeVirial,
111                tensor                     virialScaled,
112                const PbcAiuc              pbcAiuc);
113
114     /*! \brief
115      * Update data-structures (e.g. after NB search step).
116      *
117      * Updates the constraints data and copies it to the GPU. Should be
118      * called if the particles were sorted, redistributed between domains, etc.
119      * Does not recycle the data preparation routines from the CPU version.
120      * All three atoms from single water molecule should be handled by the same GPU.
121      *
122      * SETTLEs atom ID's is taken from idef.il[F_SETTLE].iatoms.
123      *
124      * \param[in] idef    System topology
125      */
126     void set(const InteractionDefinitions& idef);
127
128 private:
129     //! GPU context object
130     const DeviceContext& deviceContext_;
131     //! GPU stream
132     const DeviceStream& deviceStream_;
133
134     //! Scaled virial tensor (9 reals, GPU)
135     std::vector<float> h_virialScaled_;
136     //! Scaled virial tensor (9 reals, GPU)
137     float* d_virialScaled_;
138
139     //! Number of settles
140     int numSettles_ = 0;
141
142     //! Indexes of atoms (.x for oxygen, .y and.z for hydrogens, CPU)
143     std::vector<int3> h_atomIds_;
144     //! Indexes of atoms (.x for oxygen, .y and.z for hydrogens, GPU)
145     int3* d_atomIds_;
146     //! Current size of the array of atom IDs
147     int numAtomIds_ = -1;
148     //! Allocated size for the array of atom IDs
149     int numAtomIdsAlloc_ = -1;
150
151     //! Settle parameters
152     SettleParameters settleParameters_;
153 };
154
155 } // namespace gmx
156
157 #endif // GMX_MDLIB_SETTLE_GPU_CUH