2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2019, 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 CUDA implementation of SETTLE
39 * \author Artem Zhmurov <zhmurov@gmail.com>
41 * \ingroup module_mdlib
43 #ifndef GMX_MDLIB_SETTLE_CUDA_CUH
44 #define GMX_MDLIB_SETTLE_CUDA_CUH
48 #include "gromacs/gpu_utils/gputraits.cuh"
49 #include "gromacs/math/functions.h"
50 #include "gromacs/math/invertmatrix.h"
51 #include "gromacs/math/vec.h"
52 #include "gromacs/mdtypes/mdatom.h"
53 #include "gromacs/pbcutil/pbc.h"
54 #include "gromacs/pbcutil/pbc_aiuc.h"
55 #include "gromacs/topology/idef.h"
56 #include "gromacs/topology/topology.h"
61 /* \brief Parameters for SETTLE algorithm.
63 * Following structure and subroutine are copy-pasted from the CPU version of SETTLE.
64 * This is a temporary solution and they will be recycled in more appropriate way.
65 * \todo Remove duplicates, check if recomputing makes more sense in some cases.
66 * \todo Move the projection parameters into separate structure.
68 struct SettleParameters
70 //! Mass of oxygen atom
72 //! Mass of hydrogen atom
74 //! Relative hydrogen mass (i.e. mH/(mO+2*mH))
76 //! Target distance between oxygen and hydrogen
78 //! Target distance between hydrogen atoms
80 //! Double relative H mass times height of the H-O-H triangle.
82 //! Height of the H-O-H triangle minus ra.
84 //! Half the H-H distance.
86 //! Reciprocal H-H distance
89 /* Parameters below are used for projection */
90 //! Reciprocal oxygen mass
92 //! Reciprocal hydrogen mass
94 //! Reciprocal O-H distance
96 //! Reciprocal H-H distance (again)
98 //! Reciprocal projection matrix
102 /*! \brief Initializes a projection matrix.
104 * \todo This is identical to to CPU code. Unification is needed.
106 * \param[in] invmO Reciprocal oxygen mass
107 * \param[in] invmH Reciprocal hydrogen mass
108 * \param[in] dOH Target O-H bond length
109 * \param[in] dHH Target H-H bond length
110 * \param[out] inverseCouplingMatrix Inverse bond coupling matrix for the projection version of SETTLE
112 static void initializeProjectionMatrix(const real invmO,
116 matrix inverseCouplingMatrix)
118 // We normalize the inverse masses with invmO for the matrix inversion.
119 // so we can keep using masses of almost zero for frozen particles,
120 // without running out of the float range in invertMatrix.
121 double invmORelative = 1.0;
122 double invmHRelative = invmH / static_cast<double>(invmO);
123 double distanceRatio = dHH / static_cast<double>(dOH);
125 /* Construct the constraint coupling matrix */
127 mat[0][0] = invmORelative + invmHRelative;
128 mat[0][1] = invmORelative * (1.0 - 0.5 * gmx::square(distanceRatio));
129 mat[0][2] = invmHRelative * 0.5 * distanceRatio;
130 mat[1][1] = mat[0][0];
131 mat[1][2] = mat[0][2];
132 mat[2][2] = invmHRelative + invmHRelative;
133 mat[1][0] = mat[0][1];
134 mat[2][0] = mat[0][2];
135 mat[2][1] = mat[1][2];
137 gmx::invertMatrix(mat, inverseCouplingMatrix);
139 msmul(inverseCouplingMatrix, 1 / invmO, inverseCouplingMatrix);
142 /*! \brief Initializes settle parameters.
144 * \todo This is identical to the CPU code. Unification is needed.
146 * \param[out] p SettleParameters structure to initialize
147 * \param[in] mO Mass of oxygen atom
148 * \param[in] mH Mass of hydrogen atom
149 * \param[in] dOH Target O-H bond length
150 * \param[in] dHH Target H-H bond length
152 gmx_unused // Temporary solution to keep clang happy
154 initSettleParameters(SettleParameters* p, const real mO, const real mH, const real dOH, const real dHH)
156 // We calculate parameters in double precision to minimize errors.
157 // The velocity correction applied during SETTLE coordinate constraining
158 // introduces a systematic error of approximately 1 bit per atom,
159 // depending on what the compiler does with the code.
162 real invmO = 1.0 / mO;
163 real invmH = 1.0 / mH;
167 wohh = mO + 2.0 * mH;
171 double rc = dHH / 2.0;
172 double ra = 2.0 * mH * std::sqrt(dOH * dOH - rc * rc) / wohh;
173 p->rb = std::sqrt(dOH * dOH - rc * rc) - ra;
178 // For projection: inverse masses and coupling matrix inversion
182 p->invdOH = 1.0 / dOH;
183 p->invdHH = 1.0 / dHH;
185 initializeProjectionMatrix(invmO, invmH, dOH, dHH, p->invmat);
188 /*! \internal \brief Class with interfaces and data for CUDA version of SETTLE. */
193 /*! \brief Create SETTLE object
195 * Extracts masses for oxygen and hydrogen as well as the O-H and H-H target distances
196 * from the topology data (mtop), check their values for consistency and calls the
197 * following constructor.
199 * \param[in] mtop Topology of the system to gen the masses for O and H atoms and
200 * target O-H and H-H distances. These values are also checked for
202 * \param[in] commandStream Device stream to use.
204 SettleCuda(const gmx_mtop_t& mtop, CommandStream commandStream);
208 /*! \brief Apply SETTLE.
210 * Applies SETTLE to coordinates and velocities, stored on GPU. Data at pointers d_xp and
211 * d_v change in the GPU memory. The results are not automatically copied back to the CPU
212 * memory. Method uses this class data structures which should be updated when needed using
215 * \param[in] d_x Coordinates before timestep (in GPU memory)
216 * \param[in,out] d_xp Coordinates after timestep (in GPU memory). The
217 * resulting constrained coordinates will be saved here.
218 * \param[in] updateVelocities If the velocities should be updated.
219 * \param[in,out] d_v Velocities to update (in GPU memory, can be nullptr
221 * \param[in] invdt Reciprocal timestep (to scale Lagrange
222 * multipliers when velocities are updated)
223 * \param[in] computeVirial If virial should be updated.
224 * \param[in,out] virialScaled Scaled virial tensor to be updated.
226 void apply(const float3* d_x,
228 const bool updateVelocities,
231 const bool computeVirial,
232 tensor virialScaled);
235 * Update data-structures (e.g. after NB search step).
237 * Updates the constraints data and copies it to the GPU. Should be
238 * called if the particles were sorted, redistributed between domains, etc.
239 * Does not recycle the data preparation routines from the CPU version.
240 * All three atoms from single water molecule should be handled by the same GPU.
242 * SETTLEs atom ID's is taken from idef.il[F_SETTLE].iatoms.
244 * \param[in] idef System topology
245 * \param[in] md Atoms data. Can be used to update masses if needed (not used now).
247 void set(const t_idef& idef, const t_mdatoms& md);
252 * Converts pbc data from t_pbc into the PbcAiuc format and stores the latter.
254 * \todo PBC should not be handled by constraints.
256 * \param[in] pbc The PBC data in t_pbc format.
258 void setPbc(const t_pbc* pbc);
263 CommandStream commandStream_;
264 //! Periodic boundary data
267 //! Scaled virial tensor (9 reals, GPU)
268 std::vector<float> h_virialScaled_;
269 //! Scaled virial tensor (9 reals, GPU)
270 float* d_virialScaled_;
272 //! Number of settles
275 //! Indexes of atoms (.x for oxygen, .y and.z for hydrogens, CPU)
276 std::vector<int3> h_atomIds_;
277 //! Indexes of atoms (.x for oxygen, .y and.z for hydrogens, GPU)
279 //! Current size of the array of atom IDs
280 int numAtomIds_ = -1;
281 //! Allocated size for the array of atom IDs
282 int numAtomIdsAlloc_ = -1;
284 //! Settle parameters
285 SettleParameters settleParameters_;