StatePropagatorDataGpu object to manage GPU forces, positions and velocities buffers
[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/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"
56
57 namespace gmx
58 {
59
60 /*! \internal \brief Class with interfaces and data for CUDA version of Update-Constraint. */
61 class UpdateConstrainCuda::Impl
62 {
63
64     public:
65         /*! \brief Create Update-Constrain object
66          *
67          * \param[in] ir             Input record data: LINCS takes number of iterations and order of
68          *                           projection from it.
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.
72          */
73         Impl(const t_inputrec     &ir,
74              const gmx_mtop_t     &mtop,
75              const void           *commandStream);
76
77         ~Impl();
78
79         /*! \brief Integrate
80          *
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.
89          *
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
99          */
100         void integrate(real                              dt,
101                        bool                              updateVelocities,
102                        bool                              computeVirial,
103                        tensor                            virial,
104                        bool                              doTempCouple,
105                        gmx::ArrayRef<const t_grp_tcstat> tcstat,
106                        bool                              doPressureCouple,
107                        float                             dtPressureCouple,
108                        const matrix                      velocityScalingMatrix);
109
110         /*! \brief Set the pointers and update data-structures (e.g. after NB search step).
111          *
112          * \param[in,out]  d_x            Device buffer with coordinates.
113          * \param[in,out]  d_v            Device buffer with velocities.
114          * \param[in]      d_f            Device buffer with forces.
115          * \param[in] idef                System topology
116          * \param[in] md                  Atoms data.
117          * \param[in] numTempScaleValues  Number of temperature scaling groups. Set zero for no temperature coupling.
118          */
119         void set(DeviceBuffer<float>        d_x,
120                  DeviceBuffer<float>        d_v,
121                  const DeviceBuffer<float>  d_f,
122                  const t_idef              &idef,
123                  const t_mdatoms           &md,
124                  const int                  numTempScaleValues);
125
126         /*! \brief
127          * Update PBC data.
128          *
129          * Converts PBC data from t_pbc into the PbcAiuc format and stores the latter.
130          *
131          * \param[in] pbc The PBC data in t_pbc format.
132          */
133         void setPbc(const t_pbc *pbc);
134
135     private:
136
137         //! CUDA stream
138         CommandStream       commandStream_         = nullptr;
139
140         //! Periodic boundary data
141         PbcAiuc             pbcAiuc_;
142
143         //! Number of atoms
144         int                 numAtoms_;
145
146         //! Local copy of the pointer to the device positions buffer
147         float3             *d_x_;
148         //! Local copy of the pointer to the device velocities buffer
149         float3             *d_v_;
150         //! Local copy of the pointer to the device forces buffer
151         float3             *d_f_;
152
153         //! Device buffer for intermediate positions (maintained internally)
154         float3             *d_xp_;
155         //! Number of elements in shifted coordinates buffer
156         int                 numXp_                 = -1;
157         //! Allocation size for the shifted coordinates buffer
158         int                 numXpAlloc_            = -1;
159
160
161         //! 1/mass for all atoms (GPU)
162         real               *d_inverseMasses_;
163         //! Number of elements in reciprocal masses buffer
164         int                 numInverseMasses_      = -1;
165         //! Allocation size for the reciprocal masses buffer
166         int                 numInverseMassesAlloc_ = -1;
167
168         //! Leap-Frog integrator
169         std::unique_ptr<LeapFrogCuda>        integrator_;
170         //! LINCS CUDA object to use for non-water constraints
171         std::unique_ptr<LincsCuda>           lincsCuda_;
172         //! SETTLE CUDA object for water constrains
173         std::unique_ptr<SettleCuda>          settleCuda_;
174
175 };
176
177 } // namespace gmx
178
179 #endif