Use ListOfLists for atom to constraint mapping
[gromacs.git] / src / gromacs / pbcutil / pbc_aiuc_cuda.cuh
blob0444a9e3e6ec129e2d3ca42f8c278ddee9cdaf95
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 Basic routines to handle periodic boundary conditions with CUDA.
38  *
39  * This file contains GPU implementation of the PBC-aware vector evaluation.
40  *
41  * \todo CPU, GPU and SIMD routines essentially do the same operations on
42  *       different data-types. Currently this leads to code duplication,
43  *       which has to be resolved. For details, see Redmine task #2863
44  *       https://redmine.gromacs.org/issues/2863
45  *
46  * \author Mark Abraham <mark.j.abraham@gmail.com>
47  * \author Berk Hess <hess@kth.se>
48  * \author Artem Zhmurov <zhmurov@gmail.com>
49  *
50  * \inlibraryapi
51  * \ingroup module_pbcutil
52  */
53 #ifndef GMX_PBCUTIL_PBC_AIUC_CUDA_CUH
54 #define GMX_PBCUTIL_PBC_AIUC_CUDA_CUH
56 #include "gromacs/gpu_utils/gpu_vec.cuh"
57 #include "gromacs/gpu_utils/vectype_ops.cuh"
58 #include "gromacs/pbcutil/pbc_aiuc.h"
60 /*! \brief Computes the vector between two points taking PBC into account.
61  *
62  * Computes the vector dr between points r2 and r1, taking into account the
63  * periodic boundary conditions, described in pbcAiuc object. Note that this
64  * routine always does the PBC arithmetic for all directions, multiplying the
65  * displacements by zeroes if the corresponding direction is not periodic.
66  * For triclinic boxes only distances up to half the smallest box diagonal
67  * element are guaranteed to be the shortest. This means that distances from
68  * 0.5/sqrt(2) times a box vector length (e.g. for a rhombic dodecahedron)
69  * can use a more distant periodic image.
70  *
71  * \todo This routine uses CUDA float4 types for input coordinates and
72  *       returns in rvec data-type. Other than that, it does essentially
73  *       the same thing as the version below, as well as SIMD and CPU
74  *       versions. This routine is used in gpubonded module.
75  *       To avoid code duplication, these implementations should be
76  *       unified. See Redmine task #2863:
77  *       https://redmine.gromacs.org/issues/2863
78  *
79  * \param[in]  pbcAiuc  PBC object.
80  * \param[in]  r1       Coordinates of the first point.
81  * \param[in]  r2       Coordinates of the second point.
82  * \param[out] dr       Resulting distance.
83  */
84 template<bool returnShift>
85 static __forceinline__ __device__ int
86                        pbcDxAiuc(const PbcAiuc& pbcAiuc, const float4& r1, const float4& r2, fvec dr)
88     dr[XX] = r1.x - r2.x;
89     dr[YY] = r1.y - r2.y;
90     dr[ZZ] = r1.z - r2.z;
92     float shz = rintf(dr[ZZ] * pbcAiuc.invBoxDiagZ);
93     dr[XX] -= shz * pbcAiuc.boxZX;
94     dr[YY] -= shz * pbcAiuc.boxZY;
95     dr[ZZ] -= shz * pbcAiuc.boxZZ;
97     float shy = rintf(dr[YY] * pbcAiuc.invBoxDiagY);
98     dr[XX] -= shy * pbcAiuc.boxYX;
99     dr[YY] -= shy * pbcAiuc.boxYY;
101     float shx = rintf(dr[XX] * pbcAiuc.invBoxDiagX);
102     dr[XX] -= shx * pbcAiuc.boxXX;
104     if (returnShift)
105     {
106         ivec ishift;
108         ishift[XX] = -__float2int_rn(shx);
109         ishift[YY] = -__float2int_rn(shy);
110         ishift[ZZ] = -__float2int_rn(shz);
112         return IVEC2IS(ishift);
113     }
114     else
115     {
116         return 0;
117     }
120 /*! \brief Computes the vector between two points taking PBC into account.
122  * Computes the vector dr between points r2 and r1, taking into account the
123  * periodic boundary conditions, described in pbcAiuc object. Same as above,
124  * only takes and returns data in float3 format. Does not return shifts.
126  * \todo This routine uses CUDA float3 types for both input and returns
127  *       values. Other than that, it does essentially the same thing as the
128  *       version above, as well as SIMD and CPU versions. This routine is
129  *       used in GPU-based constraints.
130  *       To avoid code duplication, these implementations should be
131  *       unified. See Redmine task #2863:
132  *       https://redmine.gromacs.org/issues/2863
134  * \param[in]  pbcAiuc  PBC object.
135  * \param[in]  r1       Coordinates of the first point.
136  * \param[in]  r2       Coordinates of the second point.
137  * \returns    dr       Resulting distance.
138  */
139 static __forceinline__ __host__ __device__ float3 pbcDxAiuc(const PbcAiuc& pbcAiuc,
140                                                             const float3&  r1,
141                                                             const float3&  r2)
143     float3 dr = r1 - r2;
145     float shz = rintf(dr.z * pbcAiuc.invBoxDiagZ);
146     dr.x -= shz * pbcAiuc.boxZX;
147     dr.y -= shz * pbcAiuc.boxZY;
148     dr.z -= shz * pbcAiuc.boxZZ;
150     float shy = rintf(dr.y * pbcAiuc.invBoxDiagY);
151     dr.x -= shy * pbcAiuc.boxYX;
152     dr.y -= shy * pbcAiuc.boxYY;
154     float shx = rintf(dr.x * pbcAiuc.invBoxDiagX);
155     dr.x -= shx * pbcAiuc.boxXX;
157     return dr;
160 #endif // GMX_PBCUTIL_PBC_AIUC_CUDA_CUH