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