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 CUDA implementation class for update and constraints.
39 * This header file is needed to include from both the device-side
40 * kernels file, and the host-side management code.
42 * \author Artem Zhmurov <zhmurov@gmail.com>
44 * \ingroup module_mdlib
46 #ifndef GMX_MDLIB_UPDATE_CONSTRAIN_CUDA_IMPL_H
47 #define GMX_MDLIB_UPDATE_CONSTRAIN_CUDA_IMPL_H
51 #include "gromacs/mdlib/leapfrog_cuda.cuh"
52 #include "gromacs/mdlib/lincs_cuda.cuh"
53 #include "gromacs/mdlib/settle_cuda.cuh"
54 #include "gromacs/mdlib/update_constrain_cuda.h"
55 #include "gromacs/mdtypes/inputrec.h"
60 /*! \internal \brief Class with interfaces and data for CUDA version of Update-Constraint. */
61 class UpdateConstrainCuda::Impl
65 /*! \brief Create Update-Constrain object
67 * \param[in] ir Input record data: LINCS takes number of iterations and order of
69 * \param[in] mtop Topology of the system: SETTLE gets the masses for O and H atoms
70 * and target O-H and H-H distances from this object.
71 * \param[in] commandStream GPU stream to use. Can be nullptr.
73 Impl(const t_inputrec &ir,
74 const gmx_mtop_t &mtop,
75 const void *commandStream);
81 * Integrates the equation of motion using Leap-Frog algorithm and applies
82 * LINCS and SETTLE constraints.
83 * Updates d_xp_ and d_v_ fields of this object.
84 * If computeVirial is true, constraints virial is written at the provided pointer.
85 * doTempCouple should be true if:
86 * 1. The temperature coupling is enabled.
87 * 2. This is the temperature coupling step.
88 * Parameters virial/lambdas can be nullptr if computeVirial/doTempCouple are false.
90 * \param[in] dt Timestep
91 * \param[in] updateVelocities If the velocities should be constrained.
92 * \param[in] computeVirial If virial should be updated.
93 * \param[out] virial Place to save virial tensor.
94 * \param[in] doTempCouple If the temperature coupling should be performed.
95 * \param[in] tcstat Temperature coupling data.
96 * \param[in] doPressureCouple If the temperature coupling should be applied.
97 * \param[in] dtPressureCouple Period between pressure coupling steps
98 * \param[in] velocityScalingMatrix Parrinello-Rahman velocity scaling matrix
100 void integrate(const real dt,
101 const bool updateVelocities,
102 const bool computeVirial,
104 const bool doTempCouple,
105 gmx::ArrayRef<const t_grp_tcstat> tcstat,
106 const bool doPressureCouple,
107 const float dtPressureCouple,
108 const matrix velocityScalingMatrix);
111 * Update data-structures (e.g. after NB search step).
113 * \param[in] idef System topology
114 * \param[in] md Atoms data.
115 * \param[in] numTempScaleValues Number of temperature scaling groups. Set zero for no temperature coupling.
117 void set(const t_idef &idef,
119 const int numTempScaleValues);
124 * Converts PBC data from t_pbc into the PbcAiuc format and stores the latter.
126 * \param[in] pbc The PBC data in t_pbc format.
128 void setPbc(const t_pbc *pbc);
131 * Copy coordinates from CPU to GPU.
133 * The data are assumed to be in float3/fvec format (single precision).
135 * \param[in] h_x CPU pointer where coordinates should be copied from.
137 void copyCoordinatesToGpu(const rvec *h_x);
140 * Copy velocities from CPU to GPU.
142 * The data are assumed to be in float3/fvec format (single precision).
144 * \param[in] h_v CPU pointer where velocities should be copied from.
146 void copyVelocitiesToGpu(const rvec *h_v);
149 * Copy forces from CPU to GPU.
151 * The data are assumed to be in float3/fvec format (single precision).
153 * \param[in] h_f CPU pointer where forces should be copied from.
155 void copyForcesToGpu(const rvec *h_f);
158 * Copy coordinates from GPU to CPU.
160 * The data are assumed to be in float3/fvec format (single precision).
162 * \param[out] h_xp CPU pointer where coordinates should be copied to.
164 void copyCoordinatesFromGpu(rvec *h_xp);
167 * Copy velocities from GPU to CPU.
169 * The velocities are assumed to be in float3/fvec format (single precision).
171 * \param[in] h_v Pointer to velocities data.
173 void copyVelocitiesFromGpu(rvec *h_v);
176 * Copy forces from GPU to CPU.
178 * The forces are assumed to be in float3/fvec format (single precision).
180 * \param[in] h_f Pointer to forces data.
182 void copyForcesFromGpu(rvec *h_f);
185 * Set the internal GPU-memory x, xprime and v pointers.
187 * Data is not copied. The data are assumed to be in float3/fvec format
188 * (float3 is used internally, but the data layout should be identical).
190 * \param[in] d_x Pointer to the coordinates for the input (on GPU)
191 * \param[in] d_xp Pointer to the coordinates for the output (on GPU)
192 * \param[in] d_v Pointer to the velocities (on GPU)
193 * \param[in] d_f Pointer to the forces (on GPU)
195 void setXVFPointers(rvec *d_x, rvec *d_xp, rvec *d_v, rvec *d_f);
200 CommandStream commandStream_ = nullptr;
202 //! Periodic boundary data
208 //! Coordinates before the timestep (on GPU)
210 //! Number of elements in coordinates buffer
212 //! Allocation size for the coordinates buffer
215 //! Coordinates after the timestep (on GPU).
217 //! Number of elements in shifted coordinates buffer
219 //! Allocation size for the shifted coordinates buffer
220 int numXpAlloc_ = -1;
222 //! Velocities of atoms (on GPU)
224 //! Number of elements in velocities buffer
226 //! Allocation size for the velocities buffer
229 //! Forces, exerted by atoms (on GPU)
231 //! Number of elements in forces buffer
233 //! Allocation size for the forces buffer
236 //! 1/mass for all atoms (GPU)
237 real *d_inverseMasses_;
238 //! Number of elements in reciprocal masses buffer
239 int numInverseMasses_ = -1;
240 //! Allocation size for the reciprocal masses buffer
241 int numInverseMassesAlloc_ = -1;
243 //! Leap-Frog integrator
244 std::unique_ptr<LeapFrogCuda> integrator_;
245 //! LINCS CUDA object to use for non-water constraints
246 std::unique_ptr<LincsCuda> lincsCuda_;
247 //! SETTLE CUDA object for water constrains
248 std::unique_ptr<SettleCuda> settleCuda_;