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