Centralize management of coordinates ready on GPU event
[alexxy/gromacs.git] / src / gromacs / mdlib / update_constrain_cuda_impl.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 /*! \internal \file
36  *
37  * \brief Declares CUDA implementation class for update and constraints.
38  *
39  * This header file is needed to include from both the device-side
40  * kernels file, and the host-side management code.
41  *
42  * \author Artem Zhmurov <zhmurov@gmail.com>
43  *
44  * \ingroup module_mdlib
45  */
46 #ifndef GMX_MDLIB_UPDATE_CONSTRAIN_CUDA_IMPL_H
47 #define GMX_MDLIB_UPDATE_CONSTRAIN_CUDA_IMPL_H
48
49 #include "gmxpre.h"
50
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"
57
58 namespace gmx
59 {
60
61 /*! \internal \brief Class with interfaces and data for CUDA version of Update-Constraint. */
62 class UpdateConstrainCuda::Impl
63 {
64
65     public:
66         /*! \brief Create Update-Constrain object.
67          *
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.
73          *
74          * \param[in] ir                Input record data: LINCS takes number of iterations and order of
75          *                              projection from it.
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.
80          */
81         Impl(const t_inputrec     &ir,
82              const gmx_mtop_t     &mtop,
83              const void           *commandStream,
84              GpuEventSynchronizer *xUpdatedOnDevice);
85
86         ~Impl();
87
88         /*! \brief Integrate
89          *
90          * Integrates the equation of motion using Leap-Frog algorithm and applies
91          * LINCS and SETTLE constraints.
92          * Updates d_xp_ and d_v_ fields of this object.
93          * If computeVirial is true, constraints virial is written at the provided pointer.
94          * doTempCouple should be true if:
95          *   1. The temperature coupling is enabled.
96          *   2. This is the temperature coupling step.
97          * Parameters virial/lambdas can be nullptr if computeVirial/doTempCouple are false.
98          *
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
108          */
109         void integrate(real                              dt,
110                        bool                              updateVelocities,
111                        bool                              computeVirial,
112                        tensor                            virial,
113                        bool                              doTempCouple,
114                        gmx::ArrayRef<const t_grp_tcstat> tcstat,
115                        bool                              doPressureCouple,
116                        float                             dtPressureCouple,
117                        const matrix                      velocityScalingMatrix);
118
119         /*! \brief Set the pointers and update data-structures (e.g. after NB search step).
120          *
121          * \param[in,out]  d_x            Device buffer with coordinates.
122          * \param[in,out]  d_v            Device buffer with velocities.
123          * \param[in]      d_f            Device buffer with forces.
124          * \param[in] idef                System topology
125          * \param[in] md                  Atoms data.
126          * \param[in] numTempScaleValues  Number of temperature scaling groups. Set zero for no temperature coupling.
127          */
128         void set(DeviceBuffer<float>        d_x,
129                  DeviceBuffer<float>        d_v,
130                  const DeviceBuffer<float>  d_f,
131                  const t_idef              &idef,
132                  const t_mdatoms           &md,
133                  const int                  numTempScaleValues);
134
135         /*! \brief
136          * Update PBC data.
137          *
138          * Converts PBC data from t_pbc into the PbcAiuc format and stores the latter.
139          *
140          * \param[in] pbc The PBC data in t_pbc format.
141          */
142         void setPbc(const t_pbc *pbc);
143
144         /*! \brief Blocking wait on the update of coordinates being ready.
145          *
146          * \todo Remove when the "stitching" is done.
147          */
148         void waitCoordinatesReadyOnDevice();
149
150
151         /*! \brief Return the synchronizer associated with the event indicated that the coordinates are ready on the device.
152          */
153         GpuEventSynchronizer* getCoordinatesReadySync();
154
155     private:
156
157         //! CUDA stream
158         CommandStream       commandStream_         = nullptr;
159
160         //! Periodic boundary data
161         PbcAiuc             pbcAiuc_;
162
163         //! Number of atoms
164         int                 numAtoms_;
165
166         //! Local copy of the pointer to the device positions buffer
167         float3             *d_x_;
168         //! Local copy of the pointer to the device velocities buffer
169         float3             *d_v_;
170         //! Local copy of the pointer to the device forces buffer
171         float3             *d_f_;
172
173         //! Device buffer for intermediate positions (maintained internally)
174         float3             *d_xp_;
175         //! Number of elements in shifted coordinates buffer
176         int                 numXp_                 = -1;
177         //! Allocation size for the shifted coordinates buffer
178         int                 numXpAlloc_            = -1;
179
180
181         //! 1/mass for all atoms (GPU)
182         real               *d_inverseMasses_;
183         //! Number of elements in reciprocal masses buffer
184         int                 numInverseMasses_      = -1;
185         //! Allocation size for the reciprocal masses buffer
186         int                 numInverseMassesAlloc_ = -1;
187
188         //! Leap-Frog integrator
189         std::unique_ptr<LeapFrogCuda>        integrator_;
190         //! LINCS CUDA object to use for non-water constraints
191         std::unique_ptr<LincsCuda>           lincsCuda_;
192         //! SETTLE CUDA object for water constrains
193         std::unique_ptr<SettleCuda>          settleCuda_;
194
195         //! An pointer to the event to indicate when the update of coordinates is complete
196         GpuEventSynchronizer                *coordinatesReady_;
197 };
198
199 } // namespace gmx
200
201 #endif