2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2013,2014, 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>
43 * \ingroup module_simd
46 #ifndef GMX_SIMD_VECTOR_OPERATIONS_H
47 #define GMX_SIMD_VECTOR_OPERATIONS_H
49 #include "gromacs/simd/simd.h"
52 /*! \addtogroup module_simd */
55 #if (defined GMX_SIMD_HAVE_FLOAT) || (defined DOXYGEN)
56 /*! \brief SIMD float inner product of multiple float vectors.
58 * For normal usage you should always call the real-precision \ref gmx_simd_iprod_r.
60 * \param ax X components of first vectors
61 * \param ay Y components of first vectors
62 * \param az Z components of first vectors
63 * \param bx X components of second vectors
64 * \param by Y components of second vectors
65 * \param bz Z components of second vectors
67 * \return Element i will be res[i] = ax[i]*bx[i]+ay[i]*by[i]+az[i]*bz[i].
69 * \note The SIMD part is that we calculate many scalar products in one call.
71 static gmx_inline gmx_simd_float_t gmx_simdcall
72 gmx_simd_iprod_f(gmx_simd_float_t ax, gmx_simd_float_t ay, gmx_simd_float_t az,
73 gmx_simd_float_t bx, gmx_simd_float_t by, gmx_simd_float_t bz)
77 ret = gmx_simd_mul_f(ax, bx);
78 ret = gmx_simd_fmadd_f(ay, by, ret);
79 ret = gmx_simd_fmadd_f(az, bz, ret);
84 /*! \brief SIMD float norm squared of multiple vectors.
86 * For normal usage you should always call the real-precision \ref gmx_simd_norm2_r.
88 * \param ax X components of vectors
89 * \param ay Y components of vectors
90 * \param az Z components of vectors
92 * \return Element i will be res[i] = ax[i]*ax[i]+ay[i]*ay[i]+az[i]*az[i].
94 * \note This corresponds to the scalar product of the vector with itself, but
95 * the compiler might be able to optimize it better with identical vectors.
97 static gmx_inline gmx_simd_float_t gmx_simdcall
98 gmx_simd_norm2_f(gmx_simd_float_t ax, gmx_simd_float_t ay, gmx_simd_float_t az)
100 gmx_simd_float_t ret;
102 ret = gmx_simd_mul_f(ax, ax);
103 ret = gmx_simd_fmadd_f(ay, ay, ret);
104 ret = gmx_simd_fmadd_f(az, az, ret);
109 /*! \brief Calculating r^2 is the same as evaluating the norm of dx*dx.
111 * For details, see \ref gmx_simd_norm2_f.
113 #define gmx_simd_calc_rsq_f gmx_simd_norm2_f
115 /*! \brief SIMD float cross-product of multiple vectors.
117 * For normal usage you should always call the real-precision \ref gmx_simd_cprod_r.
119 * \param ax X components of first vectors
120 * \param ay Y components of first vectors
121 * \param az Z components of first vectors
122 * \param bx X components of second vectors
123 * \param by Y components of second vectors
124 * \param bz Z components of second vectors
125 * \param[out] cx X components of cross product vectors
126 * \param[out] cy Y components of cross product vectors
127 * \param[out] cz Z components of cross product vectors
131 * This calculates C = A x B, where the cross denotes the cross product.
132 * The arguments x/y/z denotes the different components, and each element
133 * corresponds to a separate vector.
135 static gmx_inline void gmx_simdcall
136 gmx_simd_cprod_f(gmx_simd_float_t ax, gmx_simd_float_t ay, gmx_simd_float_t az,
137 gmx_simd_float_t bx, gmx_simd_float_t by, gmx_simd_float_t bz,
138 gmx_simd_float_t *cx, gmx_simd_float_t *cy, gmx_simd_float_t *cz)
140 *cx = gmx_simd_mul_f(ay, bz);
141 *cx = gmx_simd_fnmadd_f(az, by, *cx);
143 *cy = gmx_simd_mul_f(az, bx);
144 *cy = gmx_simd_fnmadd_f(ax, bz, *cy);
146 *cz = gmx_simd_mul_f(ax, by);
147 *cz = gmx_simd_fnmadd_f(ay, bx, *cz);
149 #endif /* GMX_SIMD_HAVE_FLOAT */
151 #if (defined GMX_SIMD_HAVE_DOUBLE) || (defined DOXYGEN)
152 /*! \brief SIMD double inner product of multiple double vectors.
154 * \copydetails gmx_simd_iprod_f
156 static gmx_inline gmx_simd_double_t gmx_simdcall
157 gmx_simd_iprod_d(gmx_simd_double_t ax, gmx_simd_double_t ay, gmx_simd_double_t az,
158 gmx_simd_double_t bx, gmx_simd_double_t by, gmx_simd_double_t bz)
160 gmx_simd_double_t ret;
162 ret = gmx_simd_mul_d(ax, bx);
163 ret = gmx_simd_fmadd_d(ay, by, ret);
164 ret = gmx_simd_fmadd_d(az, bz, ret);
169 /*! \brief SIMD double norm squared of multiple vectors.
171 * \copydetails gmx_simd_norm2_f
173 static gmx_inline gmx_simd_double_t gmx_simdcall
174 gmx_simd_norm2_d(gmx_simd_double_t ax, gmx_simd_double_t ay, gmx_simd_double_t az)
176 gmx_simd_double_t ret;
178 ret = gmx_simd_mul_d(ax, ax);
179 ret = gmx_simd_fmadd_d(ay, ay, ret);
180 ret = gmx_simd_fmadd_d(az, az, ret);
185 /*! \brief Calculating r^2 is the same as evaluating the norm of dx*dx.
187 * For details, see \ref gmx_simd_norm2_d.
189 #define gmx_simd_calc_rsq_d gmx_simd_norm2_d
191 /*! \brief SIMD double cross-product of multiple vectors.
193 * \copydetails gmx_simd_cprod_f
195 static gmx_inline void gmx_simdcall
196 gmx_simd_cprod_d(gmx_simd_double_t ax, gmx_simd_double_t ay, gmx_simd_double_t az,
197 gmx_simd_double_t bx, gmx_simd_double_t by, gmx_simd_double_t bz,
198 gmx_simd_double_t *cx, gmx_simd_double_t *cy, gmx_simd_double_t *cz)
200 *cx = gmx_simd_mul_d(ay, bz);
201 *cx = gmx_simd_fnmadd_d(az, by, *cx);
203 *cy = gmx_simd_mul_d(az, bx);
204 *cy = gmx_simd_fnmadd_d(ax, bz, *cy);
206 *cz = gmx_simd_mul_d(ax, by);
207 *cz = gmx_simd_fnmadd_d(ay, bx, *cz);
209 #endif /* GMX_SIMD_HAVE_DOUBLE */
212 #if (defined GMX_SIMD4_HAVE_FLOAT) || (defined DOXYGEN)
213 /*! \brief SIMD4 float inner product of four float vectors.
215 * \copydetails gmx_simd_norm2_f
217 static gmx_inline gmx_simd4_float_t gmx_simdcall
218 gmx_simd4_norm2_f(gmx_simd4_float_t ax, gmx_simd4_float_t ay, gmx_simd4_float_t az)
220 gmx_simd4_float_t ret;
222 ret = gmx_simd4_mul_f(ax, ax);
223 ret = gmx_simd4_fmadd_f(ay, ay, ret);
224 ret = gmx_simd4_fmadd_f(az, az, ret);
229 /*! \brief Calculating r^2 is the same as evaluating the norm of dx*dx.
231 * For details, see \ref gmx_simd4_norm2_f
233 #define gmx_simd4_calc_rsq_f gmx_simd4_norm2_f
235 #endif /* GMX_SIMD4_HAVE_FLOAT */
237 #if (defined GMX_SIMD4_HAVE_DOUBLE) || (defined DOXYGEN)
238 /*! \brief SIMD4 double norm squared of multiple vectors.
240 * \copydetails gmx_simd_norm2_f
242 static gmx_inline gmx_simd4_double_t gmx_simdcall
243 gmx_simd4_norm2_d(gmx_simd4_double_t ax, gmx_simd4_double_t ay, gmx_simd4_double_t az)
245 gmx_simd4_double_t ret;
247 ret = gmx_simd4_mul_d(ax, ax);
248 ret = gmx_simd4_fmadd_d(ay, ay, ret);
249 ret = gmx_simd4_fmadd_d(az, az, ret);
254 /*! \brief Calculating r^2 is the same as evaluating the norm of dx*dx.
256 * For details, see \ref gmx_simd4_norm2_d.
258 #define gmx_simd4_calc_rsq_d gmx_simd4_norm2_d
260 #endif /* GMX_SIMD4_HAVE_DOUBLE */
264 /* Documented for the single branch below */
265 # define gmx_simd_iprod_r gmx_simd_iprod_d
266 # define gmx_simd_norm2_r gmx_simd_norm2_d
267 # define gmx_simd_calc_rsq_r gmx_simd_calc_rsq_d
268 # define gmx_simd_cprod_r gmx_simd_cprod_d
269 # define gmx_simd4_norm2_r gmx_simd4_norm2_d
270 # define gmx_simd4_calc_rsq_r gmx_simd4_calc_rsq_d
271 #else /* GMX_DOUBLE */
273 /*! \brief SIMD real inner product of multiple real vectors.
275 * This will call \ref gmx_simd_iprod_d if GMX_DOUBLE is defined, otherwise
276 * \ref gmx_simd_iprod_f.
278 * \copydetails gmx_simd_iprod_f
280 # define gmx_simd_iprod_r gmx_simd_iprod_f
282 /*! \brief SIMD real norm squared of multiple real vectors.
284 * This will call \ref gmx_simd_norm2_d if GMX_DOUBLE is defined, otherwise
285 * \ref gmx_simd_norm2_f.
287 * \copydetails gmx_simd_norm2_f
289 # define gmx_simd_norm2_r gmx_simd_norm2_f
291 /*! \brief Calculating r^2 is the same as evaluating the norm of dx*dx.
293 * This will call \ref gmx_simd_calc_rsq_d if GMX_DOUBLE is defined, otherwise
294 * \ref gmx_simd_calc_rsq_f.
296 * \copydetails gmx_simd_calc_rsq_f
298 # define gmx_simd_calc_rsq_r gmx_simd_calc_rsq_f
300 /*! \brief SIMD real cross-product of multiple real vectors.
302 * This will call \ref gmx_simd_cprod_d if GMX_DOUBLE is defined, otherwise
303 * \ref gmx_simd_cprod_f.
305 * \copydetails gmx_simd_cprod_f
307 # define gmx_simd_cprod_r gmx_simd_cprod_f
309 /*! \brief SIMD4 real norm squared of multiple vectors.
311 * This will call \ref gmx_simd4_norm2_d if GMX_DOUBLE is defined, otherwise
312 * \ref gmx_simd4_norm2_f.
314 * \copydetails gmx_simd4_norm2_f
316 # define gmx_simd4_norm2_r gmx_simd4_norm2_f
318 /*! \brief Calculating r^2 is the same as evaluating the norm of dx*dx.
320 * This will call \ref gmx_simd4_calc_rsq_d if GMX_DOUBLE is defined, otherwise
321 * \ref gmx_simd4_calc_rsq_f.
323 * \copydetails gmx_simd4_calc_rsq_f
325 # define gmx_simd4_calc_rsq_r gmx_simd4_calc_rsq_f
327 #endif /* GMX_DOUBLE */
332 #endif /* GMX_SIMD_VECTOR_OPERATIONS_H */