2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2014,2015,2016,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.
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.
35 #ifndef GMX_SIMD_IMPL_ARM_NEON_SIMD_FLOAT_H
36 #define GMX_SIMD_IMPL_ARM_NEON_SIMD_FLOAT_H
46 #include "gromacs/math/utilities.h"
56 SimdFloat(float f) : simdInternal_(vdupq_n_f32(f)) {}
58 // Internal utility constructor to simplify return statements
59 SimdFloat(float32x4_t simd) : simdInternal_(simd) {}
61 float32x4_t simdInternal_;
69 SimdFInt32(std::int32_t i) : simdInternal_(vdupq_n_s32(i)) {}
71 // Internal utility constructor to simplify return statements
72 SimdFInt32(int32x4_t simd) : simdInternal_(simd) {}
74 int32x4_t simdInternal_;
82 SimdFBool(bool b) : simdInternal_(vdupq_n_u32( b ? 0xFFFFFFFF : 0)) {}
84 // Internal utility constructor to simplify return statements
85 SimdFBool(uint32x4_t simd) : simdInternal_(simd) {}
87 uint32x4_t simdInternal_;
95 SimdFIBool(bool b) : simdInternal_(vdupq_n_u32( b ? 0xFFFFFFFF : 0)) {}
97 // Internal utility constructor to simplify return statements
98 SimdFIBool(uint32x4_t simd) : simdInternal_(simd) {}
100 uint32x4_t simdInternal_;
103 static inline SimdFloat gmx_simdcall
104 simdLoad(const float *m, SimdFloatTag = {})
106 assert(std::size_t(m) % 16 == 0);
112 static inline void gmx_simdcall
113 store(float *m, SimdFloat a)
115 assert(std::size_t(m) % 16 == 0);
116 vst1q_f32(m, a.simdInternal_);
119 static inline SimdFloat gmx_simdcall
120 simdLoadU(const float *m, SimdFloatTag = {})
127 static inline void gmx_simdcall
128 storeU(float *m, SimdFloat a)
130 vst1q_f32(m, a.simdInternal_);
133 static inline SimdFloat gmx_simdcall
141 static inline SimdFInt32 gmx_simdcall
142 simdLoad(const std::int32_t * m, SimdFInt32Tag)
144 assert(std::size_t(m) % 16 == 0);
150 static inline void gmx_simdcall
151 store(std::int32_t * m, SimdFInt32 a)
153 assert(std::size_t(m) % 16 == 0);
154 vst1q_s32(m, a.simdInternal_);
157 static inline SimdFInt32 gmx_simdcall
158 simdLoadU(const std::int32_t *m, SimdFInt32Tag)
165 static inline void gmx_simdcall
166 storeU(std::int32_t * m, SimdFInt32 a)
168 vst1q_s32(m, a.simdInternal_);
171 static inline SimdFInt32 gmx_simdcall
179 template<int index> gmx_simdcall
180 static inline std::int32_t
181 extract(SimdFInt32 a)
183 return vgetq_lane_s32(a.simdInternal_, index);
186 static inline SimdFloat gmx_simdcall
187 operator&(SimdFloat a, SimdFloat b)
190 vreinterpretq_f32_s32(vandq_s32(vreinterpretq_s32_f32(a.simdInternal_),
191 vreinterpretq_s32_f32(b.simdInternal_)))
195 static inline SimdFloat gmx_simdcall
196 andNot(SimdFloat a, SimdFloat b)
199 vreinterpretq_f32_s32(vbicq_s32(vreinterpretq_s32_f32(b.simdInternal_),
200 vreinterpretq_s32_f32(a.simdInternal_)))
204 static inline SimdFloat gmx_simdcall
205 operator|(SimdFloat a, SimdFloat b)
208 vreinterpretq_f32_s32(vorrq_s32(vreinterpretq_s32_f32(a.simdInternal_),
209 vreinterpretq_s32_f32(b.simdInternal_)))
213 static inline SimdFloat gmx_simdcall
214 operator^(SimdFloat a, SimdFloat b)
217 vreinterpretq_f32_s32(veorq_s32(vreinterpretq_s32_f32(a.simdInternal_),
218 vreinterpretq_s32_f32(b.simdInternal_)))
222 static inline SimdFloat gmx_simdcall
223 operator+(SimdFloat a, SimdFloat b)
226 vaddq_f32(a.simdInternal_, b.simdInternal_)
230 static inline SimdFloat gmx_simdcall
231 operator-(SimdFloat a, SimdFloat b)
234 vsubq_f32(a.simdInternal_, b.simdInternal_)
238 static inline SimdFloat gmx_simdcall
239 operator-(SimdFloat x)
242 vnegq_f32(x.simdInternal_)
246 static inline SimdFloat gmx_simdcall
247 operator*(SimdFloat a, SimdFloat b)
250 vmulq_f32(a.simdInternal_, b.simdInternal_)
254 // Override for Neon-Asimd
255 #if GMX_SIMD_ARM_NEON
256 static inline SimdFloat gmx_simdcall
257 fma(SimdFloat a, SimdFloat b, SimdFloat c)
260 #ifdef __ARM_FEATURE_FMA
261 vfmaq_f32(c.simdInternal_, b.simdInternal_, a.simdInternal_)
263 vmlaq_f32(c.simdInternal_, b.simdInternal_, a.simdInternal_)
268 static inline SimdFloat gmx_simdcall
269 fms(SimdFloat a, SimdFloat b, SimdFloat c)
272 #ifdef __ARM_FEATURE_FMA
273 vnegq_f32(vfmsq_f32(c.simdInternal_, b.simdInternal_, a.simdInternal_))
275 vnegq_f32(vmlsq_f32(c.simdInternal_, b.simdInternal_, a.simdInternal_))
280 static inline SimdFloat gmx_simdcall
281 fnma(SimdFloat a, SimdFloat b, SimdFloat c)
284 #ifdef __ARM_FEATURE_FMA
285 vfmsq_f32(c.simdInternal_, b.simdInternal_, a.simdInternal_)
287 vmlsq_f32(c.simdInternal_, b.simdInternal_, a.simdInternal_)
292 static inline SimdFloat gmx_simdcall
293 fnms(SimdFloat a, SimdFloat b, SimdFloat c)
296 #ifdef __ARM_FEATURE_FMA
297 vnegq_f32(vfmaq_f32(c.simdInternal_, b.simdInternal_, a.simdInternal_))
299 vnegq_f32(vmlaq_f32(c.simdInternal_, b.simdInternal_, a.simdInternal_))
305 static inline SimdFloat gmx_simdcall
309 vrsqrteq_f32(x.simdInternal_)
313 // The SIMD implementation seems to overflow when we square lu for
314 // values close to FLOAT_MAX, so we fall back on the version in
315 // simd_math.h, which is probably slightly slower.
316 #if GMX_SIMD_HAVE_NATIVE_RSQRT_ITER_FLOAT
317 static inline SimdFloat gmx_simdcall
318 rsqrtIter(SimdFloat lu, SimdFloat x)
321 vmulq_f32(lu.simdInternal_, vrsqrtsq_f32(vmulq_f32(lu.simdInternal_, lu.simdInternal_), x.simdInternal_))
326 static inline SimdFloat gmx_simdcall
330 vrecpeq_f32(x.simdInternal_)
334 static inline SimdFloat gmx_simdcall
335 rcpIter(SimdFloat lu, SimdFloat x)
338 vmulq_f32(lu.simdInternal_, vrecpsq_f32(lu.simdInternal_, x.simdInternal_))
342 static inline SimdFloat gmx_simdcall
343 maskAdd(SimdFloat a, SimdFloat b, SimdFBool m)
345 b.simdInternal_ = vreinterpretq_f32_u32(vandq_u32(vreinterpretq_u32_f32(b.simdInternal_),
349 vaddq_f32(a.simdInternal_, b.simdInternal_)
353 static inline SimdFloat gmx_simdcall
354 maskzMul(SimdFloat a, SimdFloat b, SimdFBool m)
356 SimdFloat tmp = a * b;
359 vreinterpretq_f32_u32(vandq_u32(vreinterpretq_u32_f32(tmp.simdInternal_),
364 static inline SimdFloat gmx_simdcall
365 maskzFma(SimdFloat a, SimdFloat b, SimdFloat c, SimdFBool m)
367 #ifdef __ARM_FEATURE_FMA
368 float32x4_t tmp = vfmaq_f32(c.simdInternal_, b.simdInternal_, a.simdInternal_);
370 float32x4_t tmp = vmlaq_f32(c.simdInternal_, b.simdInternal_, a.simdInternal_);
374 vreinterpretq_f32_u32(vandq_u32(vreinterpretq_u32_f32(tmp),
379 static inline SimdFloat gmx_simdcall
380 maskzRsqrt(SimdFloat x, SimdFBool m)
382 // The result will always be correct since we mask the result with m, but
383 // for debug builds we also want to make sure not to generate FP exceptions
385 x.simdInternal_ = vbslq_f32(m.simdInternal_, x.simdInternal_, vdupq_n_f32(1.0F));
388 vreinterpretq_f32_u32(vandq_u32(vreinterpretq_u32_f32(vrsqrteq_f32(x.simdInternal_)),
393 static inline SimdFloat gmx_simdcall
394 maskzRcp(SimdFloat x, SimdFBool m)
396 // The result will always be correct since we mask the result with m, but
397 // for debug builds we also want to make sure not to generate FP exceptions
399 x.simdInternal_ = vbslq_f32(m.simdInternal_, x.simdInternal_, vdupq_n_f32(1.0F));
402 vreinterpretq_f32_u32(vandq_u32(vreinterpretq_u32_f32(vrecpeq_f32(x.simdInternal_)),
407 static inline SimdFloat gmx_simdcall
411 vabsq_f32( x.simdInternal_ )
415 static inline SimdFloat gmx_simdcall
416 max(SimdFloat a, SimdFloat b)
419 vmaxq_f32(a.simdInternal_, b.simdInternal_)
423 static inline SimdFloat gmx_simdcall
424 min(SimdFloat a, SimdFloat b)
427 vminq_f32(a.simdInternal_, b.simdInternal_)
431 // Round and trunc operations are defined at the end of this file, since they
432 // need to use float-to-integer and integer-to-float conversions.
434 static inline SimdFloat gmx_simdcall
435 frexp(SimdFloat value, SimdFInt32 * exponent)
437 const int32x4_t exponentMask = vdupq_n_s32(0x7F800000);
438 const int32x4_t mantissaMask = vdupq_n_s32(0x807FFFFF);
439 const int32x4_t exponentBias = vdupq_n_s32(126); // add 1 to make our definition identical to frexp()
440 const float32x4_t half = vdupq_n_f32(0.5F);
443 iExponent = vandq_s32(vreinterpretq_s32_f32(value.simdInternal_), exponentMask);
444 iExponent = vsubq_s32(vshrq_n_s32(iExponent, 23), exponentBias);
445 exponent->simdInternal_ = iExponent;
448 vreinterpretq_f32_s32(vorrq_s32(vandq_s32(vreinterpretq_s32_f32(value.simdInternal_),
450 vreinterpretq_s32_f32(half)))
454 template <MathOptimization opt = MathOptimization::Safe>
455 static inline SimdFloat gmx_simdcall
456 ldexp(SimdFloat value, SimdFInt32 exponent)
458 const int32x4_t exponentBias = vdupq_n_s32(127);
459 int32x4_t iExponent = vaddq_s32(exponent.simdInternal_, exponentBias);
461 if (opt == MathOptimization::Safe)
463 // Make sure biased argument is not negative
464 iExponent = vmaxq_s32(iExponent, vdupq_n_s32(0));
467 iExponent = vshlq_n_s32( iExponent, 23);
470 vmulq_f32(value.simdInternal_, vreinterpretq_f32_s32(iExponent))
474 // Override for Neon-Asimd
475 #if GMX_SIMD_ARM_NEON
476 static inline float gmx_simdcall
479 float32x4_t x = a.simdInternal_;
480 float32x4_t y = vextq_f32(x, x, 2);
483 y = vextq_f32(x, x, 1);
485 return vgetq_lane_f32(x, 0);
489 static inline SimdFBool gmx_simdcall
490 operator==(SimdFloat a, SimdFloat b)
493 vceqq_f32(a.simdInternal_, b.simdInternal_)
497 static inline SimdFBool gmx_simdcall
498 operator!=(SimdFloat a, SimdFloat b)
501 vmvnq_u32(vceqq_f32(a.simdInternal_, b.simdInternal_))
505 static inline SimdFBool gmx_simdcall
506 operator<(SimdFloat a, SimdFloat b)
509 vcltq_f32(a.simdInternal_, b.simdInternal_)
513 static inline SimdFBool gmx_simdcall
514 operator<=(SimdFloat a, SimdFloat b)
517 vcleq_f32(a.simdInternal_, b.simdInternal_)
521 static inline SimdFBool gmx_simdcall
522 testBits(SimdFloat a)
524 uint32x4_t tmp = vreinterpretq_u32_f32(a.simdInternal_);
531 static inline SimdFBool gmx_simdcall
532 operator&&(SimdFBool a, SimdFBool b)
536 vandq_u32(a.simdInternal_, b.simdInternal_)
540 static inline SimdFBool gmx_simdcall
541 operator||(SimdFBool a, SimdFBool b)
544 vorrq_u32(a.simdInternal_, b.simdInternal_)
548 // Override for Neon-Asimd
549 #if GMX_SIMD_ARM_NEON
550 static inline bool gmx_simdcall
553 uint32x4_t x = a.simdInternal_;
554 uint32x4_t y = vextq_u32(x, x, 2);
557 y = vextq_u32(x, x, 1);
559 return (vgetq_lane_u32(x, 0) != 0);
563 static inline SimdFloat gmx_simdcall
564 selectByMask(SimdFloat a, SimdFBool m)
567 vreinterpretq_f32_u32(vandq_u32(vreinterpretq_u32_f32(a.simdInternal_),
572 static inline SimdFloat gmx_simdcall
573 selectByNotMask(SimdFloat a, SimdFBool m)
576 vreinterpretq_f32_u32(vbicq_u32(vreinterpretq_u32_f32(a.simdInternal_),
581 static inline SimdFloat gmx_simdcall
582 blend(SimdFloat a, SimdFloat b, SimdFBool sel)
585 vbslq_f32(sel.simdInternal_, b.simdInternal_, a.simdInternal_)
589 static inline SimdFInt32 gmx_simdcall
590 operator&(SimdFInt32 a, SimdFInt32 b)
593 vandq_s32(a.simdInternal_, b.simdInternal_)
597 static inline SimdFInt32 gmx_simdcall
598 andNot(SimdFInt32 a, SimdFInt32 b)
601 vbicq_s32(b.simdInternal_, a.simdInternal_)
605 static inline SimdFInt32 gmx_simdcall
606 operator|(SimdFInt32 a, SimdFInt32 b)
609 vorrq_s32(a.simdInternal_, b.simdInternal_)
613 static inline SimdFInt32 gmx_simdcall
614 operator^(SimdFInt32 a, SimdFInt32 b)
617 veorq_s32(a.simdInternal_, b.simdInternal_)
621 static inline SimdFInt32 gmx_simdcall
622 operator+(SimdFInt32 a, SimdFInt32 b)
625 vaddq_s32(a.simdInternal_, b.simdInternal_)
629 static inline SimdFInt32 gmx_simdcall
630 operator-(SimdFInt32 a, SimdFInt32 b)
633 vsubq_s32(a.simdInternal_, b.simdInternal_)
637 static inline SimdFInt32 gmx_simdcall
638 operator*(SimdFInt32 a, SimdFInt32 b)
641 vmulq_s32(a.simdInternal_, b.simdInternal_)
645 static inline SimdFIBool gmx_simdcall
646 operator==(SimdFInt32 a, SimdFInt32 b)
649 vceqq_s32(a.simdInternal_, b.simdInternal_)
653 static inline SimdFIBool gmx_simdcall
654 testBits(SimdFInt32 a)
657 vtstq_s32(a.simdInternal_, a.simdInternal_)
661 static inline SimdFIBool gmx_simdcall
662 operator<(SimdFInt32 a, SimdFInt32 b)
665 vcltq_s32(a.simdInternal_, b.simdInternal_)
669 static inline SimdFIBool gmx_simdcall
670 operator&&(SimdFIBool a, SimdFIBool b)
673 vandq_u32(a.simdInternal_, b.simdInternal_)
677 static inline SimdFIBool gmx_simdcall
678 operator||(SimdFIBool a, SimdFIBool b)
681 vorrq_u32(a.simdInternal_, b.simdInternal_)
685 // Override for Neon-Asimd
686 #if GMX_SIMD_ARM_NEON
687 static inline bool gmx_simdcall
688 anyTrue(SimdFIBool a)
690 uint32x4_t x = a.simdInternal_;
691 uint32x4_t y = vextq_u32(x, x, 2);
694 y = vextq_u32(x, x, 1);
696 return (vgetq_lane_u32(x, 0) != 0);
700 static inline SimdFInt32 gmx_simdcall
701 selectByMask(SimdFInt32 a, SimdFIBool m)
704 vandq_s32(a.simdInternal_, vreinterpretq_s32_u32(m.simdInternal_))
708 static inline SimdFInt32 gmx_simdcall
709 selectByNotMask(SimdFInt32 a, SimdFIBool m)
712 vbicq_s32(a.simdInternal_, vreinterpretq_s32_u32(m.simdInternal_))
716 static inline SimdFInt32 gmx_simdcall
717 blend(SimdFInt32 a, SimdFInt32 b, SimdFIBool sel)
720 vbslq_s32(sel.simdInternal_, b.simdInternal_, a.simdInternal_)
724 // Override for Neon-Asimd
725 #if GMX_SIMD_ARM_NEON
726 static inline SimdFInt32 gmx_simdcall
729 float32x4_t signBitOfA = vreinterpretq_f32_u32(vandq_u32(vdupq_n_u32(0x80000000), vreinterpretq_u32_f32(a.simdInternal_)));
730 float32x4_t half = vdupq_n_f32(0.5F);
731 float32x4_t corr = vreinterpretq_f32_u32(vorrq_u32(vreinterpretq_u32_f32(half), vreinterpretq_u32_f32(signBitOfA)));
734 vcvtq_s32_f32(vaddq_f32(a.simdInternal_, corr))
739 static inline SimdFInt32 gmx_simdcall
743 vcvtq_s32_f32(a.simdInternal_)
747 static inline SimdFloat gmx_simdcall
751 vcvtq_f32_s32(a.simdInternal_)
755 static inline SimdFIBool gmx_simdcall
763 static inline SimdFBool gmx_simdcall
764 cvtIB2B(SimdFIBool a)
771 // Override for Neon-Asimd
772 #if GMX_SIMD_ARM_NEON
773 static inline SimdFloat gmx_simdcall
776 return cvtI2R(cvtR2I(x));
779 static inline SimdFloat gmx_simdcall
782 return cvtI2R(cvttR2I(x));
788 #endif // GMX_SIMD_IMPL_ARM_NEON_SIMD_FLOAT_H