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/gpu_utils/gpueventsynchronizer.cuh"
52 #include "gromacs/mdlib/leapfrog_cuda.cuh"
53 #include "gromacs/mdlib/lincs_cuda.cuh"
54 #include "gromacs/mdlib/settle_cuda.cuh"
55 #include "gromacs/mdlib/update_constrain_cuda.h"
56 #include "gromacs/mdtypes/inputrec.h"
61 /*! \internal \brief Class with interfaces and data for CUDA version of Update-Constraint. */
62 class UpdateConstrainCuda::Impl
66 /*! \brief Create Update-Constrain object.
68 * The constructor is given a non-nullptr \p commandStream, in which all the update and constrain
69 * routines are executed. \p xUpdatedOnDevice should mark the completion of all kernels that modify
70 * coordinates. The event is maintained outside this class and also passed to all (if any) consumers
71 * of the updated coordinates. The \p xUpdatedOnDevice also can not be a nullptr because the
72 * markEvent(...) method is called unconditionally.
74 * \param[in] ir Input record data: LINCS takes number of iterations and order of
76 * \param[in] mtop Topology of the system: SETTLE gets the masses for O and H atoms
77 * and target O-H and H-H distances from this object.
78 * \param[in] commandStream GPU stream to use. Can be nullptr.
79 * \param[in] xUpdatedOnDevice The event synchronizer to use to mark that update is done on the GPU.
81 Impl(const t_inputrec &ir,
82 const gmx_mtop_t &mtop,
83 const void *commandStream,
84 GpuEventSynchronizer *xUpdatedOnDevice);
90 * Integrates the equation of motion using Leap-Frog algorithm and applies
91 * LINCS and SETTLE constraints.
92 * If computeVirial is true, constraints virial is written at the provided pointer.
93 * doTempCouple should be true if:
94 * 1. The temperature coupling is enabled.
95 * 2. This is the temperature coupling step.
96 * Parameters virial/lambdas can be nullptr if computeVirial/doTempCouple are false.
98 * \param[in] fReadyOnDevice Event synchronizer indicating that the forces are ready in the device memory.
99 * \param[in] dt Timestep.
100 * \param[in] updateVelocities If the velocities should be constrained.
101 * \param[in] computeVirial If virial should be updated.
102 * \param[out] virial Place to save virial tensor.
103 * \param[in] doTempCouple If the temperature coupling should be performed.
104 * \param[in] tcstat Temperature coupling data.
105 * \param[in] doPressureCouple If the temperature coupling should be applied.
106 * \param[in] dtPressureCouple Period between pressure coupling steps
107 * \param[in] velocityScalingMatrix Parrinello-Rahman velocity scaling matrix
109 void integrate(GpuEventSynchronizer *fReadyOnDevice,
111 bool updateVelocities,
115 gmx::ArrayRef<const t_grp_tcstat> tcstat,
116 bool doPressureCouple,
117 float dtPressureCouple,
118 const matrix velocityScalingMatrix);
120 /*! \brief Set the pointers and update data-structures (e.g. after NB search step).
122 * \param[in,out] d_x Device buffer with coordinates.
123 * \param[in,out] d_v Device buffer with velocities.
124 * \param[in] d_f Device buffer with forces.
125 * \param[in] idef System topology
126 * \param[in] md Atoms data.
127 * \param[in] numTempScaleValues Number of temperature scaling groups. Set zero for no temperature coupling.
129 void set(DeviceBuffer<float> d_x,
130 DeviceBuffer<float> d_v,
131 const DeviceBuffer<float> d_f,
134 const int numTempScaleValues);
139 * Converts PBC data from t_pbc into the PbcAiuc format and stores the latter.
141 * \param[in] pbc The PBC data in t_pbc format.
143 void setPbc(const t_pbc *pbc);
145 /*! \brief Return the synchronizer associated with the event indicated that the coordinates are ready on the device.
147 GpuEventSynchronizer* getCoordinatesReadySync();
152 CommandStream commandStream_ = nullptr;
154 //! Periodic boundary data
160 //! Local copy of the pointer to the device positions buffer
162 //! Local copy of the pointer to the device velocities buffer
164 //! Local copy of the pointer to the device forces buffer
167 //! Device buffer for intermediate positions (maintained internally)
169 //! Number of elements in shifted coordinates buffer
171 //! Allocation size for the shifted coordinates buffer
172 int numXpAlloc_ = -1;
175 //! 1/mass for all atoms (GPU)
176 real *d_inverseMasses_;
177 //! Number of elements in reciprocal masses buffer
178 int numInverseMasses_ = -1;
179 //! Allocation size for the reciprocal masses buffer
180 int numInverseMassesAlloc_ = -1;
182 //! Leap-Frog integrator
183 std::unique_ptr<LeapFrogCuda> integrator_;
184 //! LINCS CUDA object to use for non-water constraints
185 std::unique_ptr<LincsCuda> lincsCuda_;
186 //! SETTLE CUDA object for water constrains
187 std::unique_ptr<SettleCuda> settleCuda_;
189 //! An pointer to the event to indicate when the update of coordinates is complete
190 GpuEventSynchronizer *coordinatesReady_;