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/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/topology.h"
57 class InteractionDefinitions;
62 /* \brief Parameters for SETTLE algorithm.
64 * Following structure and subroutine are copy-pasted from the CPU version of SETTLE.
65 * This is a temporary solution and they will be recycled in more appropriate way.
66 * \todo Remove duplicates, check if recomputing makes more sense in some cases.
67 * \todo Move the projection parameters into separate structure.
69 struct SettleParameters
71 //! Mass of oxygen atom
73 //! Mass of hydrogen atom
75 //! Relative hydrogen mass (i.e. mH/(mO+2*mH))
77 //! Target distance between oxygen and hydrogen
79 //! Target distance between hydrogen atoms
81 //! Double relative H mass times height of the H-O-H triangle.
83 //! Height of the H-O-H triangle minus ra.
85 //! Half the H-H distance.
87 //! Reciprocal H-H distance
90 /* Parameters below are used for projection */
91 //! Reciprocal oxygen mass
93 //! Reciprocal hydrogen mass
95 //! Reciprocal O-H distance
97 //! Reciprocal H-H distance (again)
99 //! Reciprocal projection matrix
103 /*! \brief Initializes a projection matrix.
105 * \todo This is identical to to CPU code. Unification is needed.
107 * \param[in] invmO Reciprocal oxygen mass
108 * \param[in] invmH Reciprocal hydrogen mass
109 * \param[in] dOH Target O-H bond length
110 * \param[in] dHH Target H-H bond length
111 * \param[out] inverseCouplingMatrix Inverse bond coupling matrix for the projection version of SETTLE
113 static void initializeProjectionMatrix(const real invmO,
117 matrix inverseCouplingMatrix)
119 // We normalize the inverse masses with invmO for the matrix inversion.
120 // so we can keep using masses of almost zero for frozen particles,
121 // without running out of the float range in invertMatrix.
122 double invmORelative = 1.0;
123 double invmHRelative = invmH / static_cast<double>(invmO);
124 double distanceRatio = dHH / static_cast<double>(dOH);
126 /* Construct the constraint coupling matrix */
128 mat[0][0] = invmORelative + invmHRelative;
129 mat[0][1] = invmORelative * (1.0 - 0.5 * gmx::square(distanceRatio));
130 mat[0][2] = invmHRelative * 0.5 * distanceRatio;
131 mat[1][1] = mat[0][0];
132 mat[1][2] = mat[0][2];
133 mat[2][2] = invmHRelative + invmHRelative;
134 mat[1][0] = mat[0][1];
135 mat[2][0] = mat[0][2];
136 mat[2][1] = mat[1][2];
138 gmx::invertMatrix(mat, inverseCouplingMatrix);
140 msmul(inverseCouplingMatrix, 1 / invmO, inverseCouplingMatrix);
143 /*! \brief Initializes settle parameters.
145 * \todo This is identical to the CPU code. Unification is needed.
147 * \param[out] p SettleParameters structure to initialize
148 * \param[in] mO Mass of oxygen atom
149 * \param[in] mH Mass of hydrogen atom
150 * \param[in] dOH Target O-H bond length
151 * \param[in] dHH Target H-H bond length
153 gmx_unused // Temporary solution to keep clang happy
155 initSettleParameters(SettleParameters* p, const real mO, const real mH, const real dOH, const real dHH)
157 // We calculate parameters in double precision to minimize errors.
158 // The velocity correction applied during SETTLE coordinate constraining
159 // introduces a systematic error of approximately 1 bit per atom,
160 // depending on what the compiler does with the code.
163 real invmO = 1.0 / mO;
164 real invmH = 1.0 / mH;
168 wohh = mO + 2.0 * mH;
172 double rc = dHH / 2.0;
173 double ra = 2.0 * mH * std::sqrt(dOH * dOH - rc * rc) / wohh;
174 p->rb = std::sqrt(dOH * dOH - rc * rc) - ra;
179 // For projection: inverse masses and coupling matrix inversion
183 p->invdOH = 1.0 / dOH;
184 p->invdHH = 1.0 / dHH;
186 initializeProjectionMatrix(invmO, invmH, dOH, dHH, p->invmat);
189 /*! \internal \brief Class with interfaces and data for GPU version of SETTLE. */
194 /*! \brief Create SETTLE object
196 * Extracts masses for oxygen and hydrogen as well as the O-H and H-H target distances
197 * from the topology data (mtop), check their values for consistency and calls the
198 * following constructor.
200 * \param[in] mtop Topology of the system to gen the masses for O and H atoms and
201 * target O-H and H-H distances. These values are also checked for
203 * \param[in] commandStream Device stream to use.
205 SettleGpu(const gmx_mtop_t& mtop, CommandStream commandStream);
209 /*! \brief Apply SETTLE.
211 * Applies SETTLE to coordinates and velocities, stored on GPU. Data at pointers d_xp and
212 * d_v change in the GPU memory. The results are not automatically copied back to the CPU
213 * memory. Method uses this class data structures which should be updated when needed using
216 * \param[in] d_x Coordinates before timestep (in GPU memory)
217 * \param[in,out] d_xp Coordinates after timestep (in GPU memory). The
218 * resulting constrained coordinates will be saved here.
219 * \param[in] updateVelocities If the velocities should be updated.
220 * \param[in,out] d_v Velocities to update (in GPU memory, can be nullptr
222 * \param[in] invdt Reciprocal timestep (to scale Lagrange
223 * multipliers when velocities are updated)
224 * \param[in] computeVirial If virial should be updated.
225 * \param[in,out] virialScaled Scaled virial tensor to be updated.
226 * \param[in] pbcAiuc PBC data.
228 void apply(const float3* d_x,
230 const bool updateVelocities,
233 const bool computeVirial,
235 const PbcAiuc pbcAiuc);
238 * Update data-structures (e.g. after NB search step).
240 * Updates the constraints data and copies it to the GPU. Should be
241 * called if the particles were sorted, redistributed between domains, etc.
242 * Does not recycle the data preparation routines from the CPU version.
243 * All three atoms from single water molecule should be handled by the same GPU.
245 * SETTLEs atom ID's is taken from idef.il[F_SETTLE].iatoms.
247 * \param[in] idef System topology
248 * \param[in] md Atoms data. Can be used to update masses if needed (not used now).
250 void set(const InteractionDefinitions& idef, const t_mdatoms& md);
254 CommandStream commandStream_;
256 //! Scaled virial tensor (9 reals, GPU)
257 std::vector<float> h_virialScaled_;
258 //! Scaled virial tensor (9 reals, GPU)
259 float* d_virialScaled_;
261 //! Number of settles
264 //! Indexes of atoms (.x for oxygen, .y and.z for hydrogens, CPU)
265 std::vector<int3> h_atomIds_;
266 //! Indexes of atoms (.x for oxygen, .y and.z for hydrogens, GPU)
268 //! Current size of the array of atom IDs
269 int numAtomIds_ = -1;
270 //! Allocated size for the array of atom IDs
271 int numAtomIdsAlloc_ = -1;
273 //! Settle parameters
274 SettleParameters settleParameters_;
279 #endif // GMX_MDLIB_SETTLE_GPU_CUH