Remove PImpl scaffolding from CUDA version of LINCS
[gromacs.git] / src / gromacs / mdlib / update_constrain_cuda_impl.h
blob1e6869c21b58410f24814bfc4294245f40c4b39a
1 /*
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.
35 /*! \internal \file
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
49 #include "gromacs/mdlib/lincs_cuda.cuh"
50 #include "gromacs/mdlib/update_constrain_cuda.h"
51 #include "gromacs/mdtypes/inputrec.h"
52 #include "gromacs/pbcutil/pbc.h"
53 #include "gromacs/pbcutil/pbc_aiuc_cuda.cuh"
54 #include "gromacs/topology/idef.h"
56 #include "leapfrog_cuda_impl.h"
57 #include "settle_cuda_impl.h"
59 namespace gmx
62 /*! \internal \brief Class with interfaces and data for CUDA version of Update-Constraint. */
63 class UpdateConstrainCuda::Impl
66 public:
67 /*! \brief Create Update-Constrain object
69 * \param[in] numAtoms Number of atoms.
70 * \param[in] ir Input record data: LINCS takes number of iterations and order of
71 * projection from it.
72 * \param[in] mtop Topology of the system: SETTLE gets the masses for O and H atoms
73 * and target O-H and H-H distances from this object.
75 Impl(int numAtoms,
76 const t_inputrec &ir,
77 const gmx_mtop_t &mtop);
79 ~Impl();
81 /*! \brief Integrate
83 * Integrates the equation of motion using Leap-Frog algorithm and applies
84 * LINCS and SETTLE constraints.
85 * Updates d_xp_ and d_v_ fields of this object.
87 * \param[in] dt Timestep
88 * \param[in] updateVelocities If the velocities should be constrained.
89 * \param[in] computeVirial If virial should be updated.
90 * \param[out] virial Place to save virial tensor.
92 void integrate(const real dt,
93 const bool updateVelocities,
94 const bool computeVirial,
95 tensor virial);
97 /*! \brief
98 * Update data-structures (e.g. after NB search step).
100 * \param[in] idef System topology
101 * \param[in] md Atoms data. Can be used to update masses if needed (not used now).
103 void set(const t_idef &idef,
104 const t_mdatoms &md);
106 /*! \brief
107 * Update PBC data.
109 * Converts PBC data from t_pbc into the PbcAiuc format and stores the latter.
111 * \param[in] pbc The PBC data in t_pbc format.
113 void setPbc(const t_pbc *pbc);
115 /*! \brief
116 * Copy coordinates from CPU to GPU.
118 * The data are assumed to be in float3/fvec format (single precision).
120 * \param[in] h_x CPU pointer where coordinates should be copied from.
122 void copyCoordinatesToGpu(const rvec *h_x);
124 /*! \brief
125 * Copy velocities from CPU to GPU.
127 * The data are assumed to be in float3/fvec format (single precision).
129 * \param[in] h_v CPU pointer where velocities should be copied from.
131 void copyVelocitiesToGpu(const rvec *h_v);
133 /*! \brief
134 * Copy forces from CPU to GPU.
136 * The data are assumed to be in float3/fvec format (single precision).
138 * \param[in] h_f CPU pointer where forces should be copied from.
140 void copyForcesToGpu(const rvec *h_f);
142 /*! \brief
143 * Copy coordinates from GPU to CPU.
145 * The data are assumed to be in float3/fvec format (single precision).
147 * \param[out] h_xp CPU pointer where coordinates should be copied to.
149 void copyCoordinatesFromGpu(rvec *h_xp);
151 /*! \brief
152 * Copy velocities from GPU to CPU.
154 * The velocities are assumed to be in float3/fvec format (single precision).
156 * \param[in] h_v Pointer to velocities data.
158 void copyVelocitiesFromGpu(rvec *h_v);
160 /*! \brief
161 * Copy forces from GPU to CPU.
163 * The forces are assumed to be in float3/fvec format (single precision).
165 * \param[in] h_f Pointer to forces data.
167 void copyForcesFromGpu(rvec *h_f);
169 /*! \brief
170 * Set the internal GPU-memory x, xprime and v pointers.
172 * Data is not copied. The data are assumed to be in float3/fvec format
173 * (float3 is used internally, but the data layout should be identical).
175 * \param[in] d_x Pointer to the coordinates for the input (on GPU)
176 * \param[in] d_xp Pointer to the coordinates for the output (on GPU)
177 * \param[in] d_v Pointer to the velocities (on GPU)
178 * \param[in] d_f Pointer to the forces (on GPU)
180 void setXVFPointers(rvec *d_x, rvec *d_xp, rvec *d_v, rvec *d_f);
182 private:
184 //! CUDA stream
185 cudaStream_t stream_;
187 //! Periodic boundary data
188 PbcAiuc pbcAiuc_;
190 //! Number of atoms
191 int numAtoms_;
193 //! Coordinates before the timestep (on GPU)
194 float3 *d_x_;
195 //! Coordinates after the timestep (on GPU).
196 float3 *d_xp_;
197 //! Velocities of atoms (on GPU)
198 float3 *d_v_;
199 //! Forces, exerted by atoms (on GPU)
200 float3 *d_f_;
202 //! 1/mass for all atoms (GPU)
203 real *d_inverseMasses_;
205 //! Leap-Frog integrator
206 std::unique_ptr<LeapFrogCuda::Impl> integrator_;
207 //! LINCS CUDA object to use for non-water constraints
208 std::unique_ptr<LincsCuda> lincsCuda_;
209 //! SETTLE CUDA object for water constrains
210 std::unique_ptr<SettleCuda::Impl> settleCuda_;
214 } // namespace gmx
216 #endif