Apply clang-format to source tree
[alexxy/gromacs.git] / src / gromacs / simd / vector_operations.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2013,2014,2015,2017,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.
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
36 /*! \libinternal \file
37  *
38  * \brief SIMD operations corresponding to Gromacs rvec datatypes.
39  *
40  * \author Erik Lindahl <erik.lindahl@scilifelab.se>
41  *
42  * \inlibraryapi
43  * \ingroup module_simd
44  */
45
46 #ifndef GMX_SIMD_VECTOR_OPERATIONS_H
47 #define GMX_SIMD_VECTOR_OPERATIONS_H
48
49 #include "config.h"
50
51 #include "gromacs/simd/simd.h"
52
53 namespace gmx
54 {
55
56 #if GMX_SIMD
57
58 /*! \cond libapi */
59 /*! \addtogroup module_simd */
60 /*! \{ */
61
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
68
69 #        if GMX_SIMD_HAVE_FLOAT || defined DOXYGEN
70 /*! \brief SIMD float inner product of multiple float vectors.
71  *
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
78  *
79  * \return Element i will be res[i] = ax[i]*bx[i]+ay[i]*by[i]+az[i]*bz[i].
80  *
81  * \note The SIMD part is that we calculate many scalar products in one call.
82  */
83 static inline SimdFloat gmx_simdcall
84                         iprod(SimdFloat ax, SimdFloat ay, SimdFloat az, SimdFloat bx, SimdFloat by, SimdFloat bz)
85 {
86     SimdFloat ret;
87
88     ret = ax * bx;
89     ret = ay * by + ret;
90     ret = az * bz + ret;
91
92     return ret;
93 }
94
95 /*! \brief SIMD float norm squared of multiple vectors.
96  *
97  * \param ax X components of vectors
98  * \param ay Y components of vectors
99  * \param az Z components of vectors
100  *
101  * \return Element i will be res[i] = ax[i]*ax[i]+ay[i]*ay[i]+az[i]*az[i].
102  *
103  * \note This corresponds to the scalar product of the vector with itself, but
104  * the compiler might be able to optimize it better with identical vectors.
105  */
106 static inline SimdFloat gmx_simdcall norm2(SimdFloat ax, SimdFloat ay, SimdFloat az)
107 {
108     SimdFloat ret;
109
110     ret = ax * ax;
111     ret = ay * ay + ret;
112     ret = az * az + ret;
113
114     return ret;
115 }
116
117 /*! \brief SIMD float cross-product of multiple vectors.
118  *
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
128  *
129  * \returns void
130  *
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.
134  */
135 static inline void gmx_simdcall cprod(SimdFloat  ax,
136                                       SimdFloat  ay,
137                                       SimdFloat  az,
138                                       SimdFloat  bx,
139                                       SimdFloat  by,
140                                       SimdFloat  bz,
141                                       SimdFloat* cx,
142                                       SimdFloat* cy,
143                                       SimdFloat* cz)
144 {
145     *cx = ay * bz;
146     *cx = fnma(az, by, *cx);
147
148     *cy = az * bx;
149     *cy = fnma(ax, bz, *cy);
150
151     *cz = ax * by;
152     *cz = fnma(ay, bx, *cz);
153 }
154 #        endif // GMX_SIMD_HAVE_FLOAT
155
156 #        if GMX_SIMD_HAVE_DOUBLE || defined DOXYGEN
157 /*! \brief SIMD double inner product of multiple double vectors.
158  *
159  * \param ax X components of first vectors
160  * \param ay Y components of first vectors
161  * \param az Z components of first vectors
162  * \param bx X components of second vectors
163  * \param by Y components of second vectors
164  * \param bz Z components of second vectors
165  *
166  * \return Element i will be res[i] = ax[i]*bx[i]+ay[i]*by[i]+az[i]*bz[i].
167  *
168  * \note The SIMD part is that we calculate many scalar products in one call.
169  */
170 static inline SimdDouble gmx_simdcall
171                          iprod(SimdDouble ax, SimdDouble ay, SimdDouble az, SimdDouble bx, SimdDouble by, SimdDouble bz)
172 {
173     SimdDouble ret;
174
175     ret = ax * bx;
176     ret = ay * by + ret;
177     ret = az * bz + ret;
178
179     return ret;
180 }
181
182 /*! \brief SIMD double norm squared of multiple vectors.
183  *
184  * \param ax X components of vectors
185  * \param ay Y components of vectors
186  * \param az Z components of vectors
187  *
188  * \return Element i will be res[i] = ax[i]*ax[i]+ay[i]*ay[i]+az[i]*az[i].
189  *
190  * \note This corresponds to the scalar product of the vector with itself, but
191  * the compiler might be able to optimize it better with identical vectors.
192  */
193 static inline SimdDouble gmx_simdcall norm2(SimdDouble ax, SimdDouble ay, SimdDouble az)
194 {
195     SimdDouble ret;
196
197     ret = ax * ax;
198     ret = ay * ay + ret;
199     ret = az * az + ret;
200
201     return ret;
202 }
203
204 /*! \brief SIMD double cross-product of multiple vectors.
205  *
206  * \param ax X components of first vectors
207  * \param ay Y components of first vectors
208  * \param az Z components of first vectors
209  * \param bx X components of second vectors
210  * \param by Y components of second vectors
211  * \param bz Z components of second vectors
212  * \param[out] cx X components of cross product vectors
213  * \param[out] cy Y components of cross product vectors
214  * \param[out] cz Z components of cross product vectors
215  *
216  * \returns void
217  *
218  * This calculates C = A x B, where the cross denotes the cross product.
219  * The arguments x/y/z denotes the different components, and each element
220  * corresponds to a separate vector.
221  */
222 static inline void gmx_simdcall cprod(SimdDouble  ax,
223                                       SimdDouble  ay,
224                                       SimdDouble  az,
225                                       SimdDouble  bx,
226                                       SimdDouble  by,
227                                       SimdDouble  bz,
228                                       SimdDouble* cx,
229                                       SimdDouble* cy,
230                                       SimdDouble* cz)
231 {
232     *cx = ay * bz;
233     *cx = *cx - az * by;
234
235     *cy = az * bx;
236     *cy = *cy - ax * bz;
237
238     *cz = ax * by;
239     *cz = *cz - ay * bx;
240 }
241 #        endif // GMX_SIMD_HAVE_DOUBLE
242
243
244 #        if GMX_SIMD4_HAVE_FLOAT || defined DOXYGEN
245 /*! \brief SIMD4 float norm squared of multiple vectors.
246  *
247  * \param ax X components of vectors
248  * \param ay Y components of vectors
249  * \param az Z components of vectors
250  *
251  * \return Element i will be res[i] = ax[i]*ax[i]+ay[i]*ay[i]+az[i]*az[i].
252  *
253  * \note This corresponds to the scalar product of the vector with itself, but
254  * the compiler might be able to optimize it better with identical vectors.
255  */
256 static inline Simd4Float gmx_simdcall norm2(Simd4Float ax, Simd4Float ay, Simd4Float az)
257 {
258     Simd4Float ret;
259
260     ret = ax * ax;
261     ret = ay * ay + ret;
262     ret = az * az + ret;
263
264     return ret;
265 }
266
267 #        endif // GMX_SIMD4_HAVE_FLOAT
268
269 #        if GMX_SIMD4_HAVE_DOUBLE || defined DOXYGEN
270 /*! \brief SIMD4 double norm squared of multiple vectors.
271  *
272  * \param ax X components of vectors
273  * \param ay Y components of vectors
274  * \param az Z components of vectors
275  *
276  * \return Element i will be res[i] = ax[i]*ax[i]+ay[i]*ay[i]+az[i]*az[i].
277  *
278  * \note This corresponds to the scalar product of the vector with itself, but
279  * the compiler might be able to optimize it better with identical vectors.
280  */
281 static inline Simd4Double gmx_simdcall norm2(Simd4Double ax, Simd4Double ay, Simd4Double az)
282 {
283     Simd4Double ret;
284
285     ret = ax * ax;
286     ret = ay * ay + ret;
287     ret = az * az + ret;
288
289     return ret;
290 }
291
292 #        endif // GMX_SIMD4_HAVE_DOUBLE
293
294 #    endif // GMX_SIMD_HAVE REAL || defined DOXYGEN
295
296 /*! \} */
297 /*! \endcond */
298
299 #endif // GMX_SIMD
300
301 } // namespace gmx
302
303 #endif // GMX_SIMD_VECTOR_OPERATIONS_H