Replace defines with constexpr in ishift
[alexxy/gromacs.git] / src / gromacs / pbcutil / pbc_aiuc_cuda.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 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 Issue #2863
44  *       https://gitlab.com/gromacs/gromacs/-/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
55
56 #include "gromacs/gpu_utils/vectype_ops.cuh"
57 #include "gromacs/pbcutil/ishift.h"
58 #include "gromacs/pbcutil/pbc_aiuc.h"
59
60 static inline __device__ int xyzToShiftIndex(int x, int y, int z)
61 {
62     return (gmx::detail::c_nBoxX * (gmx::detail::c_nBoxY * ((z) + gmx::c_dBoxZ) + (y) + gmx::c_dBoxY)
63             + (x) + gmx::c_dBoxX);
64 }
65
66 static inline __device__ int int3ToShiftIndex(int3 iv)
67 {
68     return (xyzToShiftIndex(iv.x, iv.y, iv.z));
69 }
70
71 /*! \brief Computes the vector between two points taking PBC into account.
72  *
73  * Computes the vector dr between points r2 and r1, taking into account the
74  * periodic boundary conditions, described in pbcAiuc object. Note that this
75  * routine always does the PBC arithmetic for all directions, multiplying the
76  * displacements by zeroes if the corresponding direction is not periodic.
77  * For triclinic boxes only distances up to half the smallest box diagonal
78  * element are guaranteed to be the shortest. This means that distances from
79  * 0.5/sqrt(2) times a box vector length (e.g. for a rhombic dodecahedron)
80  * can use a more distant periodic image.
81  *
82  * \todo This routine uses CUDA float4 types for input coordinates and
83  *       returns in rvec data-type. Other than that, it does essentially
84  *       the same thing as the version below, as well as SIMD and CPU
85  *       versions. This routine is used in gpubonded module.
86  *       To avoid code duplication, these implementations should be
87  *       unified. See Issue #2863:
88  *       https://gitlab.com/gromacs/gromacs/-/issues/2863
89  *
90  * \param[in]  pbcAiuc  PBC object.
91  * \param[in]  r1       Coordinates of the first point.
92  * \param[in]  r2       Coordinates of the second point.
93  * \param[out] dr       Resulting distance.
94  */
95 template<bool returnShift>
96 static __forceinline__ __device__ int
97                        pbcDxAiuc(const PbcAiuc& pbcAiuc, const float4 r1, const float4 r2, float3& dr)
98 {
99     dr.x = r1.x - r2.x;
100     dr.y = r1.y - r2.y;
101     dr.z = r1.z - r2.z;
102
103     float shz = rintf(dr.z * pbcAiuc.invBoxDiagZ);
104     dr.x -= shz * pbcAiuc.boxZX;
105     dr.y -= shz * pbcAiuc.boxZY;
106     dr.z -= shz * pbcAiuc.boxZZ;
107
108     float shy = rintf(dr.y * pbcAiuc.invBoxDiagY);
109     dr.x -= shy * pbcAiuc.boxYX;
110     dr.y -= shy * pbcAiuc.boxYY;
111
112     float shx = rintf(dr.x * pbcAiuc.invBoxDiagX);
113     dr.x -= shx * pbcAiuc.boxXX;
114
115     if (returnShift)
116     {
117         int3 ishift;
118
119         ishift.x = -__float2int_rn(shx);
120         ishift.y = -__float2int_rn(shy);
121         ishift.z = -__float2int_rn(shz);
122
123         return int3ToShiftIndex(ishift);
124     }
125     else
126     {
127         return 0;
128     }
129 }
130
131 /*! \brief Computes the vector between two points taking PBC into account.
132  *
133  * Computes the vector dr between points r2 and r1, taking into account the
134  * periodic boundary conditions, described in pbcAiuc object. Same as above,
135  * only takes and returns data in float3 format. Does not return shifts.
136  *
137  * \todo This routine uses CUDA float3 types for both input and returns
138  *       values. Other than that, it does essentially the same thing as the
139  *       version above, as well as SIMD and CPU versions. This routine is
140  *       used in GPU-based constraints.
141  *       To avoid code duplication, these implementations should be
142  *       unified. See Issue #2863:
143  *       https://gitlab.com/gromacs/gromacs/-/issues/2863
144  *
145  * \param[in]  pbcAiuc  PBC object.
146  * \param[in]  r1       Coordinates of the first point.
147  * \param[in]  r2       Coordinates of the second point.
148  * \returns    dr       Resulting distance.
149  */
150 static __forceinline__ __host__ __device__ float3 pbcDxAiuc(const PbcAiuc& pbcAiuc,
151                                                             const float3&  r1,
152                                                             const float3&  r2)
153 {
154     float3 dr = r1 - r2;
155
156     float shz = rintf(dr.z * pbcAiuc.invBoxDiagZ);
157     dr.x -= shz * pbcAiuc.boxZX;
158     dr.y -= shz * pbcAiuc.boxZY;
159     dr.z -= shz * pbcAiuc.boxZZ;
160
161     float shy = rintf(dr.y * pbcAiuc.invBoxDiagY);
162     dr.x -= shy * pbcAiuc.boxYX;
163     dr.y -= shy * pbcAiuc.boxYY;
164
165     float shx = rintf(dr.x * pbcAiuc.invBoxDiagX);
166     dr.x -= shx * pbcAiuc.boxXX;
167
168     return dr;
169 }
170
171 #endif // GMX_PBCUTIL_PBC_AIUC_CUDA_CUH