Use DeviceBuffer in GPU update and NBNXM code
[alexxy/gromacs.git] / src / gromacs / mdlib / lincs_gpu.cuh
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2019,2020,2021, 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 /*! \libinternal \file
36  *
37  * \brief Declares the class for GPU implementation of LINCS.
38  *
39  * \author Artem Zhmurov <zhmurov@gmail.com>
40  *
41  * \ingroup module_mdlib
42  * \inlibraryapi
43  */
44 #ifndef GMX_MDLIB_LINCS_GPU_CUH
45 #define GMX_MDLIB_LINCS_GPU_CUH
46
47 #include <memory>
48
49 #include "gromacs/gpu_utils/devicebuffer_datatype.h"
50 #include "gromacs/gpu_utils/device_context.h"
51 #include "gromacs/gpu_utils/device_stream.h"
52 #include "gromacs/gpu_utils/gputraits.cuh"
53 #include "gromacs/mdlib/constr.h"
54 #include "gromacs/pbcutil/pbc_aiuc.h"
55
56 class InteractionDefinitions;
57
58 namespace gmx
59 {
60
61 /* \brief LINCS parameters and GPU pointers
62  *
63  * This is used to accumulate all the parameters and pointers so they can be passed
64  * to the GPU as a single structure.
65  *
66  */
67 struct LincsGpuKernelParameters
68 {
69     //! Periodic boundary data
70     PbcAiuc pbcAiuc;
71     //! Order of expansion when inverting the matrix
72     int expansionOrder;
73     //! Number of iterations used to correct the projection
74     int numIterations;
75     //! 1/mass for all atoms (GPU)
76     float* d_inverseMasses;
77     //! Scaled virial tensor (6 floats: [XX, XY, XZ, YY, YZ, ZZ], GPU)
78     float* d_virialScaled;
79     /*! \brief Total number of threads.
80      *
81      *  This covers all constraints and gaps in the ends of the thread blocks
82      *  that are necessary to avoid inter-block synchronizations.
83      *  Should be a multiple of block size (the last block is filled with dummy to the end).
84      */
85     int numConstraintsThreads;
86     //! List of constrained atoms (GPU memory)
87     int2* d_constraints;
88     //! Equilibrium distances for the constraints (GPU)
89     float* d_constraintsTargetLengths;
90     //! Number of constraints, coupled with the current one (GPU)
91     int* d_coupledConstraintsCounts;
92     //! List of coupled with the current one (GPU)
93     int* d_coupledConstraintsIndices;
94     //! Elements of the coupling matrix.
95     float* d_matrixA;
96     //! Mass factors (GPU)
97     float* d_massFactors;
98 };
99
100 /*! \internal \brief Class with interfaces and data for GPU version of LINCS. */
101 class LincsGpu
102 {
103
104 public:
105     /*! \brief Constructor.
106      *
107      * \param[in] numIterations    Number of iteration for the correction of the projection.
108      * \param[in] expansionOrder   Order of the matrix inversion algorithm.
109      * \param[in] deviceContext    Device context (dummy in CUDA).
110      * \param[in] deviceStream     Device command stream.
111      */
112     LincsGpu(int                  numIterations,
113              int                  expansionOrder,
114              const DeviceContext& deviceContext,
115              const DeviceStream&  deviceStream);
116     /*! \brief Destructor.*/
117     ~LincsGpu();
118
119     /*! \brief Apply LINCS.
120      *
121      * Applies LINCS to coordinates and velocities, stored on GPU.
122      * The results are not automatically copied back to the CPU memory.
123      * Method uses this class data structures which should be updated
124      * when needed using set() method.
125      *
126      * \param[in]     d_x               Coordinates before timestep (in GPU memory)
127      * \param[in,out] d_xp              Coordinates after timestep (in GPU memory). The
128      *                                  resulting constrained coordinates will be saved here.
129      * \param[in]     updateVelocities  If the velocities should be updated.
130      * \param[in,out] d_v               Velocities to update (in GPU memory, can be nullptr
131      *                                  if not updated)
132      * \param[in]     invdt             Reciprocal timestep (to scale Lagrange
133      *                                  multipliers when velocities are updated)
134      * \param[in]     computeVirial     If virial should be updated.
135      * \param[in,out] virialScaled      Scaled virial tensor to be updated.
136      * \param[in]     pbcAiuc           PBC data.
137      */
138     void apply(const DeviceBuffer<Float3> d_x,
139                DeviceBuffer<Float3>       d_xp,
140                const bool                 updateVelocities,
141                DeviceBuffer<Float3>       d_v,
142                const real                 invdt,
143                const bool                 computeVirial,
144                tensor                     virialScaled,
145                const PbcAiuc              pbcAiuc);
146
147     /*! \brief
148      * Update data-structures (e.g. after NB search step).
149      *
150      * Updates the constraints data and copies it to the GPU. Should be
151      * called if the particles were sorted, redistributed between domains, etc.
152      * This version uses common data formats so it can be called from anywhere
153      * in the code. Does not recycle the data preparation routines from the CPU
154      * version. Works only with simple case when all the constraints in idef are
155      * are handled by a single GPU. Triangles are not handled as special case.
156      *
157      * Information about constraints is taken from:
158      *     idef.il[F_CONSTR].iatoms  --- type (T) of constraint and two atom indexes (i1, i2)
159      *     idef.iparams[T].constr.dA --- target length for constraint of type T
160      *
161      * \param[in] idef      Local topology data to get information on constraints from.
162      * \param[in] numAtoms  Number of atoms.
163      * \param[in] invmass   Inverse masses of atoms.
164      */
165     void set(const InteractionDefinitions& idef, int numAtoms, const real* invmass);
166
167     /*! \brief
168      * Returns whether the maximum number of coupled constraints is supported
169      * by the GPU LINCS code.
170      *
171      * \param[in] mtop The molecular topology
172      */
173     static bool isNumCoupledConstraintsSupported(const gmx_mtop_t& mtop);
174
175 private:
176     //! GPU context object
177     const DeviceContext& deviceContext_;
178     //! GPU stream
179     const DeviceStream& deviceStream_;
180
181     //! Parameters and pointers, passed to the GPU kernel
182     LincsGpuKernelParameters kernelParams_;
183
184     //! Scaled virial tensor (6 floats: [XX, XY, XZ, YY, YZ, ZZ])
185     std::vector<float> h_virialScaled_;
186
187     /*! \brief Maximum total number of constraints so far.
188      *
189      * If the new number of constraints is larger then previous maximum, the GPU data arrays are
190      * reallocated.
191      */
192     int numConstraintsThreadsAlloc_;
193
194     /*! \brief Maximum total number of atoms so far.
195      *
196      * If the new number of atoms is larger then previous maximum, the GPU array with masses is
197      * reallocated.
198      */
199     int numAtomsAlloc_;
200 };
201
202 } // namespace gmx
203
204 #endif // GMX_MDLIB_LINCS_GPU_CUH