Remove PImpl scaffolding from CUDA version of LINCS
[gromacs.git] / src / gromacs / mdlib / lincs_cuda.cuh
blob1fb9eada234f82cc900f843ca9f7953cb9e52b1d
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 /*! \libinternal \file
36  *
37  * \brief Declaration of high-level functions of CUDA implementation of LINCS.
38  *
39  * \author Artem Zhmurov <zhmurov@gmail.com>
40  *
41  * \ingroup module_mdlib
42  * \inlibraryapi
43  */
44 #ifndef GMX_MDLIB_LINCS_CUDA_CUH
45 #define GMX_MDLIB_LINCS_CUDA_CUH
47 #include "gromacs/mdlib/constr.h"
48 #include "gromacs/mdtypes/mdatom.h"
49 #include "gromacs/pbcutil/pbc_aiuc.h"
50 #include "gromacs/topology/idef.h"
51 #include "gromacs/utility/classhelpers.h"
53 namespace gmx
56 /* \brief LINCS parameters and GPU pointers
57  *
58  * This is used to accumulate all the parameters and pointers so they can be passed
59  * to the GPU as a single structure.
60  *
61  */
62 struct LincsCudaKernelParameters
64     //! Periodic boundary data
65     PbcAiuc             pbcAiuc;
66     //! Order of expansion when inverting the matrix
67     int                 expansionOrder;
68     //! Number of iterations used to correct the projection
69     int                 numIterations;
70     //! 1/mass for all atoms (GPU)
71     float              *d_inverseMasses;
72     //! Scaled virial tensor (6 floats: [XX, XY, XZ, YY, YZ, ZZ], GPU)
73     float              *d_virialScaled;
74     /*! \brief Total number of threads.
75      *
76      *  This covers all constraints and gaps in the ends of the thread blocks
77      *  that are necessary to avoid inter-block synchronizations.
78      *  Should be a multiple of block size (the last block is filled with dummy to the end).
79      */
80     int                 numConstraintsThreads;
81     //! List of constrained atoms (GPU memory)
82     int2               *d_constraints;
83     //! Equilibrium distances for the constraints (GPU)
84     float              *d_constraintsTargetLengths;
85     //! Number of constraints, coupled with the current one (GPU)
86     int                *d_coupledConstraintsCounts;
87     //! List of coupled with the current one (GPU)
88     int                *d_coupledConstraintsIndices;
89     //! Elements of the coupling matrix.
90     float              *d_matrixA;
91     //! Mass factors (GPU)
92     float              *d_massFactors;
95 /*! \internal \brief Class with interfaces and data for CUDA version of LINCS. */
96 class LincsCuda
99     public:
101         /*! \brief Constructor.
102          *
103          * \param[in] numIterations    Number of iteration for the correction of the projection.
104          * \param[in] expansionOrder   Order of the matrix inversion algorithm.
105          */
106         LincsCuda(int numIterations,
107                   int expansionOrder);
108         /*! \brief Destructor.*/
109         ~LincsCuda();
111         /*! \brief Apply LINCS.
112          *
113          * Applies LINCS to coordinates and velocities, stored on GPU.
114          * The results are not automatically copied back to the CPU memory.
115          * Method uses this class data structures which should be updated
116          * when needed using set() and setPbc() method.
117          *
118          * \param[in]     d_x               Coordinates before timestep (in GPU memory)
119          * \param[in,out] d_xp              Coordinates after timestep (in GPU memory). The
120          *                                  resulting constrained coordinates will be saved here.
121          * \param[in]     updateVelocities  If the velocities should be updated.
122          * \param[in,out] d_v               Velocities to update (in GPU memory, can be nullptr
123          *                                  if not updated)
124          * \param[in]     invdt             Reciprocal timestep (to scale Lagrange
125          *                                  multipliers when velocities are updated)
126          * \param[in]     computeVirial     If virial should be updated.
127          * \param[in,out] virialScaled      Scaled virial tensor to be updated.
128          */
129         void apply(const float3 *d_x,
130                    float3       *d_xp,
131                    const bool    updateVelocities,
132                    float3       *d_v,
133                    const real    invdt,
134                    const bool    computeVirial,
135                    tensor        virialScaled);
137         /*! \brief
138          * Update data-structures (e.g. after NB search step).
139          *
140          * Updates the constraints data and copies it to the GPU. Should be
141          * called if the particles were sorted, redistributed between domains, etc.
142          * This version uses common data formats so it can be called from anywhere
143          * in the code. Does not recycle the data preparation routines from the CPU
144          * version. Works only with simple case when all the constraints in idef are
145          * are handled by a single GPU. Triangles are not handled as special case.
146          *
147          * Information about constraints is taken from:
148          *     idef.il[F_CONSTR].iatoms  --- type (T) of constraint and two atom indexes (i1, i2)
149          *     idef.iparams[T].constr.dA --- target length for constraint of type T
150          * From t_mdatom, the code takes:
151          *     md.invmass  --- array of inverse square root of masses for each atom in the system.
152          *
153          * \param[in] idef  Local topology data to get information on constraints from.
154          * \param[in] md    Atoms data to get atom masses from.
155          */
156         void set(const t_idef    &idef,
157                  const t_mdatoms &md);
159         /*! \brief
160          * Update PBC data.
161          *
162          * Converts pbc data from t_pbc into the PbcAiuc format and stores the latter.
163          *
164          * \todo Remove this method. LINCS should not manage PBC.
165          *
166          * \param[in] pbc The PBC data in t_pbc format.
167          */
168         void setPbc(const t_pbc *pbc);
171     private:
173         /*! \brief CUDA stream
174          *
175          *  \todo Currently default stream (nullptr) is used. The streams should be managed outside this module
176          *        and passed here if needed.
177          */
178         cudaStream_t              stream_;
180         //! Parameters and pointers, passed to the CUDA kernel
181         LincsCudaKernelParameters kernelParams_;
183         //! Scaled virial tensor (6 floats: [XX, XY, XZ, YY, YZ, ZZ])
184         std::vector<float>        h_virialScaled_;
186         /*! \brief Maximum total number of constraints so far.
187          *
188          * If the new number of constraints is larger then previous maximum, the GPU data arrays are
189          * reallocated.
190          */
191         int                       numConstraintsThreadsAlloc_;
193         /*! \brief Maximum total number of atoms so far.
194          *
195          * If the new number of atoms is larger then previous maximum, the GPU array with masses is
196          * reallocated.
197          */
198         int                       numAtomsAlloc_;
202 }      // namespace gmx
204 #endif // GMX_MDLIB_LINCS_CUDA_CUH