Remove PImpl scaffolding from CUDA version of LINCS
[gromacs.git] / src / gromacs / pbcutil / pbc_aiuc.h
blob72811bcdf91151a9769c11dbd25153b1efccd85a
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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.
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 /*! \file
37 * \brief Structure and basic routines to handle periodic boundary conditions.
39 * This file contains CPU
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
46 * \author Mark Abraham <mark.j.abraham@gmail.com>
47 * \author Berk Hess <hess@kth.se>
48 * \author Artem Zhmurov <zhmurov@gmail.com>
50 * \inlibraryapi
51 * \ingroup module_pbcutil
53 #ifndef GMX_PBCUTIL_PBC_AIUC_H
54 #define GMX_PBCUTIL_PBC_AIUC_H
56 #include "gromacs/pbcutil/ishift.h"
58 /*! \brief Compact and ordered version of the PBC matrix.
60 * The structure contains all the dimensions of the periodic box,
61 * arranged so that the memory access pattern is more efficient.
62 * This duplicates the information, stored in PBC 'box' matrix object,
63 * but without duplicating off-diagonal members of the matrix.
64 * The structure can be set by setPbcAiuc( ... ) routine below.
66 struct PbcAiuc
68 //! 1/box[ZZ][ZZ]
69 float invBoxDiagZ;
70 //! box[ZZ][XX]
71 float boxZX;
72 //! box[ZZ][YY]
73 float boxZY;
74 //! box[ZZ][ZZ]
75 float boxZZ;
76 //! 1/box[YY][YY]
77 float invBoxDiagY;
78 //! box[YY][XX]
79 float boxYX;
80 //! box[YY][YY]
81 float boxYY;
82 //! 1/box[XX][XX]
83 float invBoxDiagX;
84 //! box[XX][XX]
85 float boxXX;
88 /*! \brief Set the PBC data-structure.
90 * \param[in] numPbcDim Number of periodic dimensions:
91 * 0 - no periodicity.
92 * 1 - periodicity along X-axis.
93 * 2 - periodicity in XY plane.
94 * 3 - periodicity along all dimensions.
95 * \param[in] box Matrix, describing the periodic cell.
96 * \param[out] pbcAiuc Pointer to PbcAiuc structure to be initialized.
99 static inline
100 void setPbcAiuc(int numPbcDim,
101 const matrix box,
102 PbcAiuc *pbcAiuc)
105 pbcAiuc->invBoxDiagZ = 0.0f;
106 pbcAiuc->boxZX = 0.0f;
107 pbcAiuc->boxZY = 0.0f;
108 pbcAiuc->boxZZ = 0.0f;
109 pbcAiuc->invBoxDiagY = 0.0f;
110 pbcAiuc->boxYX = 0.0f;
111 pbcAiuc->boxYY = 0.0f;
112 pbcAiuc->invBoxDiagX = 0.0f;
113 pbcAiuc->boxXX = 0.0f;
115 if (numPbcDim > ZZ)
117 pbcAiuc->invBoxDiagZ = 1.0f/box[ZZ][ZZ];
118 pbcAiuc->boxZX = box[ZZ][XX];
119 pbcAiuc->boxZY = box[ZZ][YY];
120 pbcAiuc->boxZZ = box[ZZ][ZZ];
122 if (numPbcDim > YY)
124 pbcAiuc->invBoxDiagY = 1.0f/box[YY][YY];
125 pbcAiuc->boxYX = box[YY][XX];
126 pbcAiuc->boxYY = box[YY][YY];
128 if (numPbcDim > XX)
130 pbcAiuc->invBoxDiagX = 1.0f/box[XX][XX];
131 pbcAiuc->boxXX = box[XX][XX];
135 /*! \brief Computes the vector between two points taking PBC into account.
137 * Computes the vector dr between points r2 and r1, taking into account the
138 * periodic boundary conditions, described in pbcAiuc object. Note that this
139 * routine always does the PBC arithmetic for all directions, multiplying the
140 * displacements by zeroes if the corresponding direction is not periodic.
141 * For triclinic boxes only distances up to half the smallest box diagonal
142 * element are guaranteed to be the shortest. This means that distances from
143 * 0.5/sqrt(2) times a box vector length (e.g. for a rhombic dodecahedron)
144 * can use a more distant periodic image.
146 * \todo This routine operates on rvec types and uses PbcAiuc to define
147 * periodic box, but essentially does the same thing as SIMD and GPU
148 * version. These will have to be unified in future to avoid code
149 * duplication. See Redmine task #2863:
150 * https://redmine.gromacs.org/issues/2863
152 * \param[in] pbcAiuc PBC object.
153 * \param[in] r1 Coordinates of the first point.
154 * \param[in] r2 Coordinates of the second point.
155 * \param[out] dr Resulting distance.
157 static inline
158 void pbcDxAiuc(const PbcAiuc &pbcAiuc,
159 const rvec &r1,
160 const rvec &r2,
161 rvec dr)
163 dr[XX] = r1[XX] - r2[XX];
164 dr[YY] = r1[YY] - r2[YY];
165 dr[ZZ] = r1[ZZ] - r2[ZZ];
167 float shz = rintf(dr[ZZ]*pbcAiuc.invBoxDiagZ);
168 dr[XX] -= shz*pbcAiuc.boxZX;
169 dr[YY] -= shz*pbcAiuc.boxZY;
170 dr[ZZ] -= shz*pbcAiuc.boxZZ;
172 float shy = rintf(dr[YY]*pbcAiuc.invBoxDiagY);
173 dr[XX] -= shy*pbcAiuc.boxYX;
174 dr[YY] -= shy*pbcAiuc.boxYY;
176 float shx = rintf(dr[XX]*pbcAiuc.invBoxDiagX);
177 dr[XX] -= shx*pbcAiuc.boxXX;
181 #endif //GMX_PBCUTIL_PBC_AIUC_H