2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2018,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.
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.
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.
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.
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.
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.
35 #ifndef GMX_PBCUTIL_GPU_PBC_CUH
36 #define GMX_PBCUTIL_GPU_PBC_CUH
38 #include "gromacs/pbcutil/ishift.h"
40 /*! \brief Compact and ordered version of the PBC matrix.
42 * The structure contains all the dimensions of the periodic box,
43 * arranged in a convenient order. This duplicates the information,
44 * stored in PBC 'box' matrix object. The structure can be set by
45 * setPbcAiuc( ... ) routine below.
62 /*! \brief Computes the vector between two points taking PBC into account.
64 * Computes the vector dr between points x2 and x1, taking into account the
65 * periodic boundary conditions, described in pbcAiuc object. Note that this
66 * routine always does the PBC arithmetic for all directions, multiplying the
67 * displacements by zeroes if the corresponding direction is not periodic.
68 * For triclinic boxes only distances up to half the smallest box diagonal
69 * element are guaranteed to be the shortest. This means that distances from
70 * 0.5/sqrt(2) times a box vector length (e.g. for a rhombic dodecahedron)
71 * can use a more distant periodic image.
73 * \param[in] pbcAiuc PBC object.
74 * \param[in] x1 Coordinates of the first point.
75 * \param[in] x2 Coordinates of the second point.
76 * \param[out] dx Resulting distance.
78 template <bool returnShift>
79 static __forceinline__ __device__
80 int pbcDxAiuc(const PbcAiuc &pbcAiuc,
89 float shz = rintf(dx[ZZ]*pbcAiuc.invBoxDiagZ);
90 dx[XX] -= shz*pbcAiuc.boxZX;
91 dx[YY] -= shz*pbcAiuc.boxZY;
92 dx[ZZ] -= shz*pbcAiuc.boxZZ;
94 float shy = rintf(dx[YY]*pbcAiuc.invBoxDiagY);
95 dx[XX] -= shy*pbcAiuc.boxYX;
96 dx[YY] -= shy*pbcAiuc.boxYY;
98 float shx = rintf(dx[XX]*pbcAiuc.invBoxDiagX);
99 dx[XX] -= shx*pbcAiuc.boxXX;
105 ishift[XX] = -__float2int_rn(shx);
106 ishift[YY] = -__float2int_rn(shy);
107 ishift[ZZ] = -__float2int_rn(shz);
109 return IVEC2IS(ishift);
117 /*! \brief Set the PBC data to use in GPU kernels.
119 * \param[in] numPbcDim Number of periodic dimensions:
120 * 0 - no periodicity.
121 * 1 - periodicity along X-axis.
122 * 2 - periodicity in XY plane.
123 * 3 - periodicity along all dimensions.
124 * \param[in] box Matrix, describing the periodic cell.
125 * \param[out] pbcAiuc Pointer to PbcAiuc structure to be initialized.
128 static void setPbcAiuc(int numPbcDim,
133 pbcAiuc->invBoxDiagZ = 0.0f;
134 pbcAiuc->boxZX = 0.0f;
135 pbcAiuc->boxZY = 0.0f;
136 pbcAiuc->boxZZ = 0.0f;
137 pbcAiuc->invBoxDiagY = 0.0f;
138 pbcAiuc->boxYX = 0.0f;
139 pbcAiuc->boxYY = 0.0f;
140 pbcAiuc->invBoxDiagX = 0.0f;
141 pbcAiuc->boxXX = 0.0f;
145 pbcAiuc->invBoxDiagZ = 1.0f/box[ZZ][ZZ];
146 pbcAiuc->boxZX = box[ZZ][XX];
147 pbcAiuc->boxZY = box[ZZ][YY];
148 pbcAiuc->boxZZ = box[ZZ][ZZ];
152 pbcAiuc->invBoxDiagY = 1.0f/box[YY][YY];
153 pbcAiuc->boxYX = box[YY][XX];
154 pbcAiuc->boxYY = box[YY][YY];
158 pbcAiuc->invBoxDiagX = 1.0f/box[XX][XX];
159 pbcAiuc->boxXX = box[XX][XX];