Apply clang-format to source tree
[alexxy/gromacs.git] / src / gromacs / simd / impl_arm_neon / impl_arm_neon_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 #ifndef GMX_SIMD_IMPL_ARM_NEON_SIMD4_FLOAT_H
36 #define GMX_SIMD_IMPL_ARM_NEON_SIMD4_FLOAT_H
37
38 #include "config.h"
39
40 #include <cassert>
41 #include <cstddef>
42
43 #include <arm_neon.h>
44
45 namespace gmx
46 {
47
48 class Simd4Float
49 {
50 public:
51     Simd4Float() {}
52
53     Simd4Float(float f) : simdInternal_(vdupq_n_f32(f)) {}
54
55     // Internal utility constructor to simplify return statements
56     Simd4Float(float32x4_t simd) : simdInternal_(simd) {}
57
58     float32x4_t simdInternal_;
59 };
60
61 class Simd4FBool
62 {
63 public:
64     Simd4FBool() {}
65
66     //! \brief Construct from scalar bool
67     Simd4FBool(bool b) : simdInternal_(vdupq_n_u32(b ? 0xFFFFFFFF : 0)) {}
68
69     // Internal utility constructor to simplify return statements
70     Simd4FBool(uint32x4_t simd) : simdInternal_(simd) {}
71
72     uint32x4_t simdInternal_;
73 };
74
75 static inline Simd4Float gmx_simdcall load4(const float* m)
76 {
77     assert(size_t(m) % 16 == 0);
78     return { vld1q_f32(m) };
79 }
80
81 static inline void gmx_simdcall store4(float* m, Simd4Float a)
82 {
83     assert(size_t(m) % 16 == 0);
84     vst1q_f32(m, a.simdInternal_);
85 }
86
87 static inline Simd4Float gmx_simdcall load4U(const float* m)
88 {
89     return { vld1q_f32(m) };
90 }
91
92 static inline void gmx_simdcall store4U(float* m, Simd4Float a)
93 {
94     vst1q_f32(m, a.simdInternal_);
95 }
96
97 static inline Simd4Float gmx_simdcall simd4SetZeroF()
98 {
99     return { vdupq_n_f32(0.0F) };
100 }
101
102 static inline Simd4Float gmx_simdcall operator&(Simd4Float a, Simd4Float b)
103 {
104     return { vreinterpretq_f32_s32(vandq_s32(vreinterpretq_s32_f32(a.simdInternal_),
105                                              vreinterpretq_s32_f32(b.simdInternal_))) };
106 }
107
108 static inline Simd4Float gmx_simdcall andNot(Simd4Float a, Simd4Float b)
109 {
110     return { vreinterpretq_f32_s32(vbicq_s32(vreinterpretq_s32_f32(b.simdInternal_),
111                                              vreinterpretq_s32_f32(a.simdInternal_))) };
112 }
113
114 static inline Simd4Float gmx_simdcall operator|(Simd4Float a, Simd4Float b)
115 {
116     return { vreinterpretq_f32_s32(vorrq_s32(vreinterpretq_s32_f32(a.simdInternal_),
117                                              vreinterpretq_s32_f32(b.simdInternal_))) };
118 }
119
120 static inline Simd4Float gmx_simdcall operator^(Simd4Float a, Simd4Float b)
121 {
122     return { vreinterpretq_f32_s32(veorq_s32(vreinterpretq_s32_f32(a.simdInternal_),
123                                              vreinterpretq_s32_f32(b.simdInternal_))) };
124 }
125
126 static inline Simd4Float gmx_simdcall operator+(Simd4Float a, Simd4Float b)
127 {
128     return { vaddq_f32(a.simdInternal_, b.simdInternal_) };
129 }
130
131 static inline Simd4Float gmx_simdcall operator-(Simd4Float a, Simd4Float b)
132 {
133     return { vsubq_f32(a.simdInternal_, b.simdInternal_) };
134 }
135
136 static inline Simd4Float gmx_simdcall operator-(Simd4Float x)
137 {
138     return { vnegq_f32(x.simdInternal_) };
139 }
140
141 static inline Simd4Float gmx_simdcall operator*(Simd4Float a, Simd4Float b)
142 {
143     return { vmulq_f32(a.simdInternal_, b.simdInternal_) };
144 }
145
146 // Override for Neon-Asimd
147 #if GMX_SIMD_ARM_NEON
148 static inline Simd4Float gmx_simdcall fma(Simd4Float a, Simd4Float b, Simd4Float c)
149 {
150     return {
151 #    ifdef __ARM_FEATURE_FMA
152         vfmaq_f32(c.simdInternal_, b.simdInternal_, a.simdInternal_)
153 #    else
154         vmlaq_f32(c.simdInternal_, b.simdInternal_, a.simdInternal_)
155 #    endif
156     };
157 }
158
159 static inline Simd4Float gmx_simdcall fms(Simd4Float a, Simd4Float b, Simd4Float c)
160 {
161     return {
162 #    ifdef __ARM_FEATURE_FMA
163         vnegq_f32(vfmsq_f32(c.simdInternal_, b.simdInternal_, a.simdInternal_))
164 #    else
165         vnegq_f32(vmlsq_f32(c.simdInternal_, b.simdInternal_, a.simdInternal_))
166 #    endif
167     };
168 }
169
170 static inline Simd4Float gmx_simdcall fnma(Simd4Float a, Simd4Float b, Simd4Float c)
171 {
172     return {
173 #    ifdef __ARM_FEATURE_FMA
174         vfmsq_f32(c.simdInternal_, b.simdInternal_, a.simdInternal_)
175 #    else
176         vmlsq_f32(c.simdInternal_, b.simdInternal_, a.simdInternal_)
177 #    endif
178     };
179 }
180
181 static inline Simd4Float gmx_simdcall fnms(Simd4Float a, Simd4Float b, Simd4Float c)
182 {
183     return {
184 #    ifdef __ARM_FEATURE_FMA
185         vnegq_f32(vfmaq_f32(c.simdInternal_, b.simdInternal_, a.simdInternal_))
186 #    else
187         vnegq_f32(vmlaq_f32(c.simdInternal_, b.simdInternal_, a.simdInternal_))
188 #    endif
189     };
190 }
191 #endif
192
193 static inline Simd4Float gmx_simdcall rsqrt(Simd4Float x)
194 {
195     return { vrsqrteq_f32(x.simdInternal_) };
196 }
197
198 static inline Simd4Float gmx_simdcall abs(Simd4Float x)
199 {
200     return { vabsq_f32(x.simdInternal_) };
201 }
202
203 static inline Simd4Float gmx_simdcall max(Simd4Float a, Simd4Float b)
204 {
205     return { vmaxq_f32(a.simdInternal_, b.simdInternal_) };
206 }
207
208 static inline Simd4Float gmx_simdcall min(Simd4Float a, Simd4Float b)
209 {
210     return { vminq_f32(a.simdInternal_, b.simdInternal_) };
211 }
212
213 // Override for Neon-Asimd
214 #if GMX_SIMD_ARM_NEON
215 static inline Simd4Float gmx_simdcall round(Simd4Float x)
216 {
217     // Convert x to nearest integer
218     float32x4_t signBitOfX = vreinterpretq_f32_u32(
219             vandq_u32(vdupq_n_u32(0x80000000), vreinterpretq_u32_f32(x.simdInternal_)));
220     float32x4_t half = vdupq_n_f32(0.5F);
221     float32x4_t corr = vreinterpretq_f32_u32(
222             vorrq_u32(vreinterpretq_u32_f32(half), vreinterpretq_u32_f32(signBitOfX)));
223
224     int32x4_t integerX = vcvtq_s32_f32(vaddq_f32(x.simdInternal_, corr));
225
226     // Convert back to float
227
228     return { vcvtq_f32_s32(integerX) };
229 }
230
231 static inline Simd4Float gmx_simdcall trunc(Simd4Float x)
232 {
233     return { vcvtq_f32_s32(vcvtq_s32_f32(x.simdInternal_)) };
234 }
235 #endif
236
237 static inline void gmx_simdcall transpose(Simd4Float* v0, Simd4Float* v1, Simd4Float* v2, Simd4Float* v3)
238 {
239     float32x4x2_t t0  = vuzpq_f32(v0->simdInternal_, v2->simdInternal_);
240     float32x4x2_t t1  = vuzpq_f32(v1->simdInternal_, v3->simdInternal_);
241     float32x4x2_t t2  = vtrnq_f32(t0.val[0], t1.val[0]);
242     float32x4x2_t t3  = vtrnq_f32(t0.val[1], t1.val[1]);
243     v0->simdInternal_ = t2.val[0];
244     v1->simdInternal_ = t3.val[0];
245     v2->simdInternal_ = t2.val[1];
246     v3->simdInternal_ = t3.val[1];
247 }
248
249 static inline Simd4FBool gmx_simdcall operator==(Simd4Float a, Simd4Float b)
250 {
251     return { vceqq_f32(a.simdInternal_, b.simdInternal_) };
252 }
253
254 static inline Simd4FBool gmx_simdcall operator!=(Simd4Float a, Simd4Float b)
255 {
256     return { vmvnq_u32(vceqq_f32(a.simdInternal_, b.simdInternal_)) };
257 }
258
259 static inline Simd4FBool gmx_simdcall operator<(Simd4Float a, Simd4Float b)
260 {
261     return { vcltq_f32(a.simdInternal_, b.simdInternal_) };
262 }
263
264 static inline Simd4FBool gmx_simdcall operator<=(Simd4Float a, Simd4Float b)
265 {
266     return { vcleq_f32(a.simdInternal_, b.simdInternal_) };
267 }
268
269 static inline Simd4FBool gmx_simdcall operator&&(Simd4FBool a, Simd4FBool b)
270 {
271     return { vandq_u32(a.simdInternal_, b.simdInternal_) };
272 }
273
274 static inline Simd4FBool gmx_simdcall operator||(Simd4FBool a, Simd4FBool b)
275 {
276     return { vorrq_u32(a.simdInternal_, b.simdInternal_) };
277 }
278
279 // Override for Neon-Asimd
280 #if GMX_SIMD_ARM_NEON
281 static inline bool gmx_simdcall anyTrue(Simd4FBool a)
282 {
283     uint32x4_t x = a.simdInternal_;
284     uint32x4_t y = vextq_u32(x, x, 2);
285
286     x = vorrq_u32(x, y);
287     y = vextq_u32(x, x, 1);
288     x = vorrq_u32(x, y);
289     return (vgetq_lane_u32(x, 0) != 0);
290 }
291 #endif
292
293 static inline Simd4Float gmx_simdcall selectByMask(Simd4Float a, Simd4FBool m)
294 {
295     return { vreinterpretq_f32_u32(vandq_u32(vreinterpretq_u32_f32(a.simdInternal_), m.simdInternal_)) };
296 }
297
298 static inline Simd4Float gmx_simdcall selectByNotMask(Simd4Float a, Simd4FBool m)
299 {
300     return { vreinterpretq_f32_u32(vbicq_u32(vreinterpretq_u32_f32(a.simdInternal_), m.simdInternal_)) };
301 }
302
303 static inline Simd4Float gmx_simdcall blend(Simd4Float a, Simd4Float b, Simd4FBool sel)
304 {
305     return { vbslq_f32(sel.simdInternal_, b.simdInternal_, a.simdInternal_) };
306 }
307
308 // Override for Neon-Asimd
309 #if GMX_SIMD_ARM_NEON
310 static inline float gmx_simdcall reduce(Simd4Float a)
311 {
312     float32x4_t x = a.simdInternal_;
313     float32x4_t y = vextq_f32(x, x, 2);
314
315     x = vaddq_f32(x, y);
316     y = vextq_f32(x, x, 1);
317     x = vaddq_f32(x, y);
318     return vgetq_lane_f32(x, 0);
319 }
320
321 static inline float gmx_simdcall dotProduct(Simd4Float a, Simd4Float b)
322 {
323     Simd4Float c;
324
325     c = a * b;
326     /* set 4th element to 0, then add all of them */
327     c.simdInternal_ = vsetq_lane_f32(0.0F, c.simdInternal_, 3);
328     return reduce(c);
329 }
330 #endif
331
332 } // namespace gmx
333
334 #endif // GMX_SIMD_IMPL_ARM_NEON_SIMD4_FLOAT_H