Apply clang-format to source tree
[alexxy/gromacs.git] / src / gromacs / simd / impl_ibm_vmx / impl_ibm_vmx_simd4_float.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2014,2015,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 #ifndef GMX_SIMD_IMPLEMENTATION_IBM_VMX_SIMD4_FLOAT_H
37 #define GMX_SIMD_IMPLEMENTATION_IBM_VMX_SIMD4_FLOAT_H
38
39 #include "config.h"
40
41 #include "gromacs/utility/basedefinitions.h"
42
43 #include "impl_ibm_vmx_definitions.h"
44
45 namespace gmx
46 {
47
48 class Simd4Float
49 {
50 public:
51     Simd4Float() {}
52
53     Simd4Float(float f)
54     {
55         __vector unsigned char perm;
56
57         simdInternal_ = vec_lde(0, const_cast<float*>(&f));
58         perm          = vec_lvsl(0, const_cast<float*>(&f));
59         simdInternal_ = vec_perm(simdInternal_, simdInternal_, perm);
60         simdInternal_ = vec_splat(simdInternal_, 0);
61     }
62
63     // Internal utility constructor to simplify return statements
64     Simd4Float(__vector float simd) : simdInternal_(simd) {}
65
66     __vector float simdInternal_;
67 };
68
69 class Simd4FBool
70 {
71 public:
72     Simd4FBool() {}
73
74     Simd4FBool(bool b)
75     {
76         simdInternal_ = reinterpret_cast<__vector vmxBool int>(vec_splat_u32(b ? 0xFFFFFFFF : 0));
77     }
78
79     // Internal utility constructor to simplify return statements
80     Simd4FBool(__vector vmxBool int simd) : simdInternal_(simd) {}
81
82     __vector vmxBool int simdInternal_;
83 };
84
85 static inline Simd4Float gmx_simdcall load4(const float* m)
86 {
87     return { vec_ld(0, const_cast<float*>(m)) };
88 }
89
90 static inline void gmx_simdcall store4(float* m, Simd4Float a)
91 {
92     vec_st(a.simdInternal_, 0, const_cast<float*>(m));
93 }
94
95 static inline Simd4Float gmx_simdcall simd4SetZeroF()
96 {
97     return { reinterpret_cast<__vector float>(vec_splat_u32(0)) };
98 }
99
100 static inline Simd4Float gmx_simdcall operator&(Simd4Float a, Simd4Float b)
101 {
102     return { vec_and(a.simdInternal_, b.simdInternal_) };
103 }
104
105 static inline Simd4Float gmx_simdcall andNot(Simd4Float a, Simd4Float b)
106 {
107     return { vec_andc(b.simdInternal_, a.simdInternal_) };
108 }
109
110 static inline Simd4Float gmx_simdcall operator|(Simd4Float a, Simd4Float b)
111 {
112     return { vec_or(a.simdInternal_, b.simdInternal_) };
113 }
114
115 static inline Simd4Float gmx_simdcall operator^(Simd4Float a, Simd4Float b)
116 {
117     return { vec_xor(a.simdInternal_, b.simdInternal_) };
118 }
119
120 static inline Simd4Float gmx_simdcall operator+(Simd4Float a, Simd4Float b)
121 {
122     return { vec_add(a.simdInternal_, b.simdInternal_) };
123 }
124
125 static inline Simd4Float gmx_simdcall operator-(Simd4Float a, Simd4Float b)
126 {
127     return { vec_sub(a.simdInternal_, b.simdInternal_) };
128 }
129
130 static inline Simd4Float gmx_simdcall operator-(Simd4Float x)
131 {
132     return { vec_xor(x.simdInternal_,
133                      reinterpret_cast<__vector float>(vec_sl(vec_splat_u32(-1), vec_splat_u32(-1)))) };
134 }
135
136 static inline Simd4Float gmx_simdcall operator*(Simd4Float a, Simd4Float b)
137 {
138     return { vec_madd(a.simdInternal_, b.simdInternal_,
139                       reinterpret_cast<__vector float>(vec_splat_u32(0))) };
140 }
141
142 static inline Simd4Float gmx_simdcall fma(Simd4Float a, Simd4Float b, Simd4Float c)
143 {
144     return { vec_madd(a.simdInternal_, b.simdInternal_, c.simdInternal_) };
145 }
146
147 static inline Simd4Float gmx_simdcall fms(Simd4Float a, Simd4Float b, Simd4Float c)
148 {
149     return { vec_madd(a.simdInternal_, b.simdInternal_, -c.simdInternal_) };
150 }
151
152 static inline Simd4Float gmx_simdcall fnma(Simd4Float a, Simd4Float b, Simd4Float c)
153 {
154     return { vec_nmsub(a.simdInternal_, b.simdInternal_, c.simdInternal_) };
155 }
156
157 static inline Simd4Float gmx_simdcall fnms(Simd4Float a, Simd4Float b, Simd4Float c)
158 {
159     return { -vec_madd(a.simdInternal_, b.simdInternal_, c.simdInternal_) };
160 }
161
162 static inline Simd4Float gmx_simdcall rsqrt(Simd4Float x)
163 {
164     return { vec_rsqrte(x.simdInternal_) };
165 }
166
167 static inline Simd4Float gmx_simdcall abs(Simd4Float x)
168 {
169     return { vec_abs(x.simdInternal_) };
170 }
171
172 static inline Simd4Float gmx_simdcall max(Simd4Float a, Simd4Float b)
173 {
174     return { vec_max(a.simdInternal_, b.simdInternal_) };
175 }
176
177 static inline Simd4Float gmx_simdcall min(Simd4Float a, Simd4Float b)
178 {
179     return { vec_min(a.simdInternal_, b.simdInternal_) };
180 }
181
182 static inline Simd4Float gmx_simdcall round(Simd4Float x)
183 {
184     return { vec_round(x.simdInternal_) };
185 }
186
187 static inline Simd4Float gmx_simdcall trunc(Simd4Float x)
188 {
189     return { vec_trunc(x.simdInternal_) };
190 }
191
192 static inline float gmx_simdcall dotProduct(Simd4Float a, Simd4Float b)
193 {
194     float res;
195
196     __vector float c = vec_madd(a.simdInternal_, b.simdInternal_,
197                                 reinterpret_cast<__vector float>(vec_splat_u32(0)));
198     // Keep only elements 0,1,2 by shifting in zero from right (xor of a vector with itself is 0)
199     c = vec_sld(c, vec_xor(a.simdInternal_, a.simdInternal_), 4);
200     // calculate sum
201     c = vec_add(c, vec_sld(c, c, 8));
202     c = vec_add(c, vec_sld(c, c, 4));
203     vec_ste(c, 0, &res);
204     return res;
205 }
206
207 static inline void gmx_simdcall transpose(Simd4Float* v0, Simd4Float* v1, Simd4Float* v2, Simd4Float* v3)
208 {
209     __vector float t0 = vec_mergeh(v0->simdInternal_, v2->simdInternal_);
210     __vector float t1 = vec_mergel(v0->simdInternal_, v2->simdInternal_);
211     __vector float t2 = vec_mergeh(v1->simdInternal_, v3->simdInternal_);
212     __vector float t3 = vec_mergel(v1->simdInternal_, v3->simdInternal_);
213     v0->simdInternal_ = vec_mergeh(t0, t2);
214     v1->simdInternal_ = vec_mergel(t0, t2);
215     v2->simdInternal_ = vec_mergeh(t1, t3);
216     v3->simdInternal_ = vec_mergel(t1, t3);
217 }
218
219 static inline Simd4FBool gmx_simdcall operator==(Simd4Float a, Simd4Float b)
220 {
221     return { vec_cmpeq(a.simdInternal_, b.simdInternal_) };
222 }
223
224 static inline Simd4FBool gmx_simdcall operator!=(Simd4Float a, Simd4Float b)
225 {
226     return { vec_or(vec_cmpgt(a.simdInternal_, b.simdInternal_),
227                     vec_cmplt(a.simdInternal_, b.simdInternal_)) };
228 }
229
230 static inline Simd4FBool gmx_simdcall operator<(Simd4Float a, Simd4Float b)
231 {
232     return { vec_cmplt(a.simdInternal_, b.simdInternal_) };
233 }
234
235 static inline Simd4FBool gmx_simdcall operator<=(Simd4Float a, Simd4Float b)
236 {
237     return { vec_cmple(a.simdInternal_, b.simdInternal_) };
238 }
239
240 static inline Simd4FBool gmx_simdcall operator&&(Simd4FBool a, Simd4FBool b)
241 {
242     return { vec_and(a.simdInternal_, b.simdInternal_) };
243 }
244
245 static inline Simd4FBool gmx_simdcall operator||(Simd4FBool a, Simd4FBool b)
246 {
247     return { vec_or(a.simdInternal_, b.simdInternal_) };
248 }
249
250 static inline bool gmx_simdcall anyTrue(Simd4FBool a)
251 {
252     return vec_any_ne(a.simdInternal_, reinterpret_cast<__vector vmxBool int>(vec_splat_u32(0)));
253 }
254
255 static inline Simd4Float gmx_simdcall selectByMask(Simd4Float a, Simd4FBool m)
256 {
257     return { vec_and(a.simdInternal_, reinterpret_cast<__vector float>(m.simdInternal_)) };
258 }
259
260 static inline Simd4Float gmx_simdcall selectByNotMask(Simd4Float a, Simd4FBool m)
261 {
262     return { vec_andc(a.simdInternal_, reinterpret_cast<__vector float>(m.simdInternal_)) };
263 }
264
265 static inline Simd4Float gmx_simdcall blend(Simd4Float a, Simd4Float b, Simd4FBool sel)
266 {
267     return { vec_sel(a.simdInternal_, b.simdInternal_, sel.simdInternal_) };
268 }
269
270 static inline float gmx_simdcall reduce(Simd4Float a)
271 {
272     __vector float c = a.simdInternal_;
273     float          res;
274
275     // calculate sum
276     c = vec_add(c, vec_sld(c, c, 8));
277     c = vec_add(c, vec_sld(c, c, 4));
278     vec_ste(c, 0, &res);
279     return res;
280 }
281
282 } // namespace gmx
283
284 #endif // GMX_SIMD_IMPLEMENTATION_IBM_VMX_SIMD4_FLOAT_H