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