2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
37 * \brief Declares class for GPU implementation of SETTLE
39 * \author Artem Zhmurov <zhmurov@gmail.com>
41 * \ingroup module_mdlib
43 #ifndef GMX_MDLIB_SETTLE_GPU_CUH
44 #define GMX_MDLIB_SETTLE_GPU_CUH
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"
58 class InteractionDefinitions;
63 /* \brief Parameters for SETTLE algorithm.
65 * Following structure and subroutine are copy-pasted from the CPU version of SETTLE.
66 * This is a temporary solution and they will be recycled in more appropriate way.
67 * \todo Remove duplicates, check if recomputing makes more sense in some cases.
68 * \todo Move the projection parameters into separate structure.
70 struct SettleParameters
72 //! Mass of oxygen atom
74 //! Mass of hydrogen atom
76 //! Relative hydrogen mass (i.e. mH/(mO+2*mH))
78 //! Target distance between oxygen and hydrogen
80 //! Target distance between hydrogen atoms
82 //! Double relative H mass times height of the H-O-H triangle.
84 //! Height of the H-O-H triangle minus ra.
86 //! Half the H-H distance.
88 //! Reciprocal H-H distance
91 /* Parameters below are used for projection */
92 //! Reciprocal oxygen mass
94 //! Reciprocal hydrogen mass
96 //! Reciprocal O-H distance
98 //! Reciprocal H-H distance (again)
100 //! Reciprocal projection matrix
104 /*! \brief Initializes a projection matrix.
106 * \todo This is identical to to CPU code. Unification is needed.
108 * \param[in] invmO Reciprocal oxygen mass
109 * \param[in] invmH Reciprocal hydrogen mass
110 * \param[in] dOH Target O-H bond length
111 * \param[in] dHH Target H-H bond length
112 * \param[out] inverseCouplingMatrix Inverse bond coupling matrix for the projection version of SETTLE
114 static void initializeProjectionMatrix(const real invmO,
118 matrix inverseCouplingMatrix)
120 // We normalize the inverse masses with invmO for the matrix inversion.
121 // so we can keep using masses of almost zero for frozen particles,
122 // without running out of the float range in invertMatrix.
123 double invmORelative = 1.0;
124 double invmHRelative = invmH / static_cast<double>(invmO);
125 double distanceRatio = dHH / static_cast<double>(dOH);
127 /* Construct the constraint coupling matrix */
129 mat[0][0] = invmORelative + invmHRelative;
130 mat[0][1] = invmORelative * (1.0 - 0.5 * gmx::square(distanceRatio));
131 mat[0][2] = invmHRelative * 0.5 * distanceRatio;
132 mat[1][1] = mat[0][0];
133 mat[1][2] = mat[0][2];
134 mat[2][2] = invmHRelative + invmHRelative;
135 mat[1][0] = mat[0][1];
136 mat[2][0] = mat[0][2];
137 mat[2][1] = mat[1][2];
139 gmx::invertMatrix(mat, inverseCouplingMatrix);
141 msmul(inverseCouplingMatrix, 1 / invmO, inverseCouplingMatrix);
144 /*! \brief Initializes settle parameters.
146 * \todo This is identical to the CPU code. Unification is needed.
148 * \param[out] p SettleParameters structure to initialize
149 * \param[in] mO Mass of oxygen atom
150 * \param[in] mH Mass of hydrogen atom
151 * \param[in] dOH Target O-H bond length
152 * \param[in] dHH Target H-H bond length
154 gmx_unused // Temporary solution to keep clang happy
156 initSettleParameters(SettleParameters* p, const real mO, const real mH, const real dOH, const real dHH)
158 // We calculate parameters in double precision to minimize errors.
159 // The velocity correction applied during SETTLE coordinate constraining
160 // introduces a systematic error of approximately 1 bit per atom,
161 // depending on what the compiler does with the code.
164 real invmO = 1.0 / mO;
165 real invmH = 1.0 / mH;
169 wohh = mO + 2.0 * mH;
173 double rc = dHH / 2.0;
174 double ra = 2.0 * mH * std::sqrt(dOH * dOH - rc * rc) / wohh;
175 p->rb = std::sqrt(dOH * dOH - rc * rc) - ra;
180 // For projection: inverse masses and coupling matrix inversion
184 p->invdOH = 1.0 / dOH;
185 p->invdHH = 1.0 / dHH;
187 initializeProjectionMatrix(invmO, invmH, dOH, dHH, p->invmat);
190 /*! \internal \brief Class with interfaces and data for GPU version of SETTLE. */
195 /*! \brief Create SETTLE object
197 * Extracts masses for oxygen and hydrogen as well as the O-H and H-H target distances
198 * from the topology data (mtop), check their values for consistency and calls the
199 * following constructor.
201 * \param[in] mtop Topology of the system to gen the masses for O and H atoms and
202 * target O-H and H-H distances. These values are also checked for
204 * \param[in] deviceContext Device context (dummy in CUDA).
205 * \param[in] deviceStream Device stream to use.
207 SettleGpu(const gmx_mtop_t& mtop, const DeviceContext& deviceContext, const DeviceStream& deviceStream);
211 /*! \brief Apply SETTLE.
213 * Applies SETTLE to coordinates and velocities, stored on GPU. Data at pointers d_xp and
214 * d_v change in the GPU memory. The results are not automatically copied back to the CPU
215 * memory. Method uses this class data structures which should be updated when needed using
218 * \param[in] d_x Coordinates before timestep (in GPU memory)
219 * \param[in,out] d_xp Coordinates after timestep (in GPU memory). The
220 * resulting constrained coordinates will be saved here.
221 * \param[in] updateVelocities If the velocities should be updated.
222 * \param[in,out] d_v Velocities to update (in GPU memory, can be nullptr
224 * \param[in] invdt Reciprocal timestep (to scale Lagrange
225 * multipliers when velocities are updated)
226 * \param[in] computeVirial If virial should be updated.
227 * \param[in,out] virialScaled Scaled virial tensor to be updated.
228 * \param[in] pbcAiuc PBC data.
230 void apply(const float3* d_x,
232 const bool updateVelocities,
235 const bool computeVirial,
237 const PbcAiuc pbcAiuc);
240 * Update data-structures (e.g. after NB search step).
242 * Updates the constraints data and copies it to the GPU. Should be
243 * called if the particles were sorted, redistributed between domains, etc.
244 * Does not recycle the data preparation routines from the CPU version.
245 * All three atoms from single water molecule should be handled by the same GPU.
247 * SETTLEs atom ID's is taken from idef.il[F_SETTLE].iatoms.
249 * \param[in] idef System topology
251 void set(const InteractionDefinitions& idef);
254 //! GPU context object
255 const DeviceContext& deviceContext_;
257 const DeviceStream& deviceStream_;
259 //! Scaled virial tensor (9 reals, GPU)
260 std::vector<float> h_virialScaled_;
261 //! Scaled virial tensor (9 reals, GPU)
262 float* d_virialScaled_;
264 //! Number of settles
267 //! Indexes of atoms (.x for oxygen, .y and.z for hydrogens, CPU)
268 std::vector<int3> h_atomIds_;
269 //! Indexes of atoms (.x for oxygen, .y and.z for hydrogens, GPU)
271 //! Current size of the array of atom IDs
272 int numAtomIds_ = -1;
273 //! Allocated size for the array of atom IDs
274 int numAtomIdsAlloc_ = -1;
276 //! Settle parameters
277 SettleParameters settleParameters_;
282 #endif // GMX_MDLIB_SETTLE_GPU_CUH