1fb9a142e1f83ea3ae895d0279377ce8ba54f789
[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) 2012,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.
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 /* The macros in this file are intended to be used for writing
37  * architecture-independent SIMD intrinsics code.
38  * To support a new architecture, adding macros here should be (nearly)
39  * all that is needed.
40  */
41
42 /* This file contains vector operation functions using SIMD intrinsics.
43  * gromacs/simd/macros.h should be included before including this file.
44  */
45
46 #ifndef GMX_SIMD_VECTOR_OPERATIONS_H
47 #define GMX_SIMD_VECTOR_OPERATIONS_H
48
49 #ifndef GMX_SIMD_MACROS_H
50 #error "gromacs/simd/macros.h was not included before including gromacs/simd/vector_operations.h"
51 #endif
52
53
54 /* x^2 + y^2 + z^2 */
55 static gmx_inline gmx_simd_real_t
56 gmx_simd_calc_rsq_r(gmx_simd_real_t x, gmx_simd_real_t y, gmx_simd_real_t z)
57 {
58     return gmx_simd_fmadd_r(z, z, gmx_simd_fmadd_r(y, y, gmx_simd_mul_r(x, x)));
59 }
60
61 /* inner-product of multiple vectors */
62 static gmx_inline gmx_simd_real_t
63 gmx_simd_iprod_r(gmx_simd_real_t ax, gmx_simd_real_t ay, gmx_simd_real_t az,
64                  gmx_simd_real_t bx, gmx_simd_real_t by, gmx_simd_real_t bz)
65 {
66     gmx_simd_real_t ret;
67
68     ret = gmx_simd_mul_r(ax, bx);
69     ret = gmx_simd_fmadd_r(ay, by, ret);
70     ret = gmx_simd_fmadd_r(az, bz, ret);
71
72     return ret;
73 }
74
75 /* norm squared of multiple vectors */
76 static gmx_inline gmx_simd_real_t
77 gmx_simd_norm2_r(gmx_simd_real_t ax, gmx_simd_real_t ay, gmx_simd_real_t az)
78 {
79     gmx_simd_real_t ret;
80
81     ret = gmx_simd_mul_r(ax, ax);
82     ret = gmx_simd_fmadd_r(ay, ay, ret);
83     ret = gmx_simd_fmadd_r(az, az, ret);
84
85     return ret;
86 }
87
88 /* cross-product of multiple vectors */
89 static gmx_inline void
90 gmx_simd_cprod_r(gmx_simd_real_t ax, gmx_simd_real_t ay, gmx_simd_real_t az,
91                  gmx_simd_real_t bx, gmx_simd_real_t by, gmx_simd_real_t bz,
92                  gmx_simd_real_t *cx, gmx_simd_real_t *cy, gmx_simd_real_t *cz)
93 {
94     *cx = gmx_simd_mul_r(ay, bz);
95     *cx = gmx_simd_fnmadd_r(az, by, *cx);
96
97     *cy = gmx_simd_mul_r(az, bx);
98     *cy = gmx_simd_fnmadd_r(ax, bz, *cy);
99
100     *cz = gmx_simd_mul_r(ax, by);
101     *cz = gmx_simd_fnmadd_r(ay, bx, *cz);
102 }
103
104 /* a + b + c + d (not really a vector operation, but where else put this?) */
105 static gmx_inline gmx_simd_real_t
106 gmx_simd_sum4_r(gmx_simd_real_t a, gmx_simd_real_t b, gmx_simd_real_t c, gmx_simd_real_t d)
107 {
108     return gmx_simd_add_r(gmx_simd_add_r(a, b), gmx_simd_add_r(c, d));
109 }
110
111
112 #endif