Introduce SimulatorBuilder
[gromacs.git] / src / gromacs / simd / vector_operations.h
blob9dca30e03a0a4315ea534971d1b152b29be3a746
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2013,2014,2015,2017, 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.
36 /*! \libinternal \file
38 * \brief SIMD operations corresponding to Gromacs rvec datatypes.
40 * \author Erik Lindahl <erik.lindahl@scilifelab.se>
42 * \inlibraryapi
43 * \ingroup module_simd
46 #ifndef GMX_SIMD_VECTOR_OPERATIONS_H
47 #define GMX_SIMD_VECTOR_OPERATIONS_H
49 #include "config.h"
51 #include "gromacs/simd/simd.h"
53 namespace gmx
56 #if GMX_SIMD
58 /*! \cond libapi */
59 /*! \addtogroup module_simd */
60 /*! \{ */
62 /* This check is not actually required, but it must be true if the
63 * code below actualy declares anything, and it makes it easy for
64 * check-source to know that this file depends on simd.h (though
65 * symbols like GMX_SIMD_HAVE_FLOAT are actually defined in its
66 * implementation headers). */
67 #if GMX_SIMD_HAVE_REAL || defined DOXYGEN
69 #if GMX_SIMD_HAVE_FLOAT || defined DOXYGEN
70 /*! \brief SIMD float inner product of multiple float vectors.
72 * \param ax X components of first vectors
73 * \param ay Y components of first vectors
74 * \param az Z components of first vectors
75 * \param bx X components of second vectors
76 * \param by Y components of second vectors
77 * \param bz Z components of second vectors
79 * \return Element i will be res[i] = ax[i]*bx[i]+ay[i]*by[i]+az[i]*bz[i].
81 * \note The SIMD part is that we calculate many scalar products in one call.
83 static inline SimdFloat gmx_simdcall
84 iprod(SimdFloat ax, SimdFloat ay, SimdFloat az,
85 SimdFloat bx, SimdFloat by, SimdFloat bz)
87 SimdFloat ret;
89 ret = ax * bx;
90 ret = ay * by + ret;
91 ret = az * bz + ret;
93 return ret;
96 /*! \brief SIMD float norm squared of multiple vectors.
98 * \param ax X components of vectors
99 * \param ay Y components of vectors
100 * \param az Z components of vectors
102 * \return Element i will be res[i] = ax[i]*ax[i]+ay[i]*ay[i]+az[i]*az[i].
104 * \note This corresponds to the scalar product of the vector with itself, but
105 * the compiler might be able to optimize it better with identical vectors.
107 static inline SimdFloat gmx_simdcall
108 norm2(SimdFloat ax, SimdFloat ay, SimdFloat az)
110 SimdFloat ret;
112 ret = ax * ax;
113 ret = ay * ay + ret;
114 ret = az * az + ret;
116 return ret;
119 /*! \brief SIMD float cross-product of multiple vectors.
121 * \param ax X components of first vectors
122 * \param ay Y components of first vectors
123 * \param az Z components of first vectors
124 * \param bx X components of second vectors
125 * \param by Y components of second vectors
126 * \param bz Z components of second vectors
127 * \param[out] cx X components of cross product vectors
128 * \param[out] cy Y components of cross product vectors
129 * \param[out] cz Z components of cross product vectors
131 * \returns void
133 * This calculates C = A x B, where the cross denotes the cross product.
134 * The arguments x/y/z denotes the different components, and each element
135 * corresponds to a separate vector.
137 static inline void gmx_simdcall
138 cprod(SimdFloat ax, SimdFloat ay, SimdFloat az,
139 SimdFloat bx, SimdFloat by, SimdFloat bz,
140 SimdFloat *cx, SimdFloat *cy, SimdFloat *cz)
142 *cx = ay * bz;
143 *cx = fnma(az, by, *cx);
145 *cy = az * bx;
146 *cy = fnma(ax, bz, *cy);
148 *cz = ax * by;
149 *cz = fnma(ay, bx, *cz);
151 #endif // GMX_SIMD_HAVE_FLOAT
153 #if GMX_SIMD_HAVE_DOUBLE || defined DOXYGEN
154 /*! \brief SIMD double inner product of multiple double vectors.
156 * \param ax X components of first vectors
157 * \param ay Y components of first vectors
158 * \param az Z components of first vectors
159 * \param bx X components of second vectors
160 * \param by Y components of second vectors
161 * \param bz Z components of second vectors
163 * \return Element i will be res[i] = ax[i]*bx[i]+ay[i]*by[i]+az[i]*bz[i].
165 * \note The SIMD part is that we calculate many scalar products in one call.
167 static inline SimdDouble gmx_simdcall
168 iprod(SimdDouble ax, SimdDouble ay, SimdDouble az,
169 SimdDouble bx, SimdDouble by, SimdDouble bz)
171 SimdDouble ret;
173 ret = ax * bx;
174 ret = ay * by + ret;
175 ret = az * bz + ret;
177 return ret;
180 /*! \brief SIMD double norm squared of multiple vectors.
182 * \param ax X components of vectors
183 * \param ay Y components of vectors
184 * \param az Z components of vectors
186 * \return Element i will be res[i] = ax[i]*ax[i]+ay[i]*ay[i]+az[i]*az[i].
188 * \note This corresponds to the scalar product of the vector with itself, but
189 * the compiler might be able to optimize it better with identical vectors.
191 static inline SimdDouble gmx_simdcall
192 norm2(SimdDouble ax, SimdDouble ay, SimdDouble az)
194 SimdDouble ret;
196 ret = ax * ax;
197 ret = ay * ay + ret;
198 ret = az * az + ret;
200 return ret;
203 /*! \brief SIMD double cross-product of multiple vectors.
205 * \param ax X components of first vectors
206 * \param ay Y components of first vectors
207 * \param az Z components of first vectors
208 * \param bx X components of second vectors
209 * \param by Y components of second vectors
210 * \param bz Z components of second vectors
211 * \param[out] cx X components of cross product vectors
212 * \param[out] cy Y components of cross product vectors
213 * \param[out] cz Z components of cross product vectors
215 * \returns void
217 * This calculates C = A x B, where the cross denotes the cross product.
218 * The arguments x/y/z denotes the different components, and each element
219 * corresponds to a separate vector.
221 static inline void gmx_simdcall
222 cprod(SimdDouble ax, SimdDouble ay, SimdDouble az,
223 SimdDouble bx, SimdDouble by, SimdDouble bz,
224 SimdDouble *cx, SimdDouble *cy, SimdDouble *cz)
226 *cx = ay * bz;
227 *cx = *cx - az * by;
229 *cy = az * bx;
230 *cy = *cy - ax * bz;
232 *cz = ax * by;
233 *cz = *cz - ay * bx;
235 #endif // GMX_SIMD_HAVE_DOUBLE
238 #if GMX_SIMD4_HAVE_FLOAT || defined DOXYGEN
239 /*! \brief SIMD4 float norm squared of multiple vectors.
241 * \param ax X components of vectors
242 * \param ay Y components of vectors
243 * \param az Z components of vectors
245 * \return Element i will be res[i] = ax[i]*ax[i]+ay[i]*ay[i]+az[i]*az[i].
247 * \note This corresponds to the scalar product of the vector with itself, but
248 * the compiler might be able to optimize it better with identical vectors.
250 static inline Simd4Float gmx_simdcall
251 norm2(Simd4Float ax, Simd4Float ay, Simd4Float az)
253 Simd4Float ret;
255 ret = ax * ax;
256 ret = ay * ay + ret;
257 ret = az * az + ret;
259 return ret;
262 #endif // GMX_SIMD4_HAVE_FLOAT
264 #if GMX_SIMD4_HAVE_DOUBLE || defined DOXYGEN
265 /*! \brief SIMD4 double norm squared of multiple vectors.
267 * \param ax X components of vectors
268 * \param ay Y components of vectors
269 * \param az Z components of vectors
271 * \return Element i will be res[i] = ax[i]*ax[i]+ay[i]*ay[i]+az[i]*az[i].
273 * \note This corresponds to the scalar product of the vector with itself, but
274 * the compiler might be able to optimize it better with identical vectors.
276 static inline Simd4Double gmx_simdcall
277 norm2(Simd4Double ax, Simd4Double ay, Simd4Double az)
279 Simd4Double ret;
281 ret = ax * ax;
282 ret = ay * ay + ret;
283 ret = az * az + ret;
285 return ret;
288 #endif // GMX_SIMD4_HAVE_DOUBLE
290 #endif // GMX_SIMD_HAVE REAL || defined DOXYGEN
292 /*! \} */
293 /*! \endcond */
295 #endif // GMX_SIMD
297 } // namespace gmx
299 #endif // GMX_SIMD_VECTOR_OPERATIONS_H