2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2014,2015,2016,2017,2019,2020, 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.
36 #ifndef GMX_SIMD_IMPLEMENTATION_IBM_VMX_SIMD_FLOAT_H
37 #define GMX_SIMD_IMPLEMENTATION_IBM_VMX_SIMD_FLOAT_H
43 #include "gromacs/math/utilities.h"
44 #include "gromacs/utility/basedefinitions.h"
46 #include "impl_ibm_vmx_definitions.h"
58 __vector
unsigned char perm
;
60 simdInternal_
= vec_lde(0, const_cast<float*>(&f
));
61 perm
= vec_lvsl(0, const_cast<float*>(&f
));
62 simdInternal_
= vec_perm(simdInternal_
, simdInternal_
, perm
);
63 simdInternal_
= vec_splat(simdInternal_
, 0);
66 // Internal utility constructor to simplify return statements
67 SimdFloat(__vector
float simd
) : simdInternal_(simd
) {}
69 __vector
float simdInternal_
;
77 SimdFInt32(std::int32_t i
)
79 __vector
unsigned char perm
;
81 simdInternal_
= vec_lde(0, const_cast<int*>(&i
));
82 perm
= vec_lvsl(0, const_cast<int*>(&i
));
83 simdInternal_
= vec_perm(simdInternal_
, simdInternal_
, perm
);
84 simdInternal_
= vec_splat(simdInternal_
, 0);
88 // Internal utility constructor to simplify return statements
89 SimdFInt32(__vector
signed int simd
) : simdInternal_(simd
) {}
91 __vector
signed int simdInternal_
;
100 simdInternal_(reinterpret_cast<__vector vmxBool
int>(vec_splats(b
? 0xFFFFFFFF : 0)))
104 // Internal utility constructor to simplify return statements
105 SimdFBool(__vector vmxBool
int simd
) : simdInternal_(simd
) {}
107 __vector vmxBool
int simdInternal_
;
116 simdInternal_(reinterpret_cast<__vector vmxBool
int>(vec_splat_u32(b
? 0xFFFFFFFF : 0)))
120 // Internal utility constructor to simplify return statements
121 SimdFIBool(__vector vmxBool
int simd
) : simdInternal_(simd
) {}
123 __vector vmxBool
int simdInternal_
;
126 static inline SimdFloat gmx_simdcall
simdLoad(const float* m
, SimdFloatTag
= {})
128 return { vec_ld(0, const_cast<float*>(m
)) };
131 static inline void gmx_simdcall
store(float* m
, SimdFloat a
)
133 vec_st(a
.simdInternal_
, 0, const_cast<float*>(m
));
136 static inline SimdFloat gmx_simdcall
setZeroF()
138 return { reinterpret_cast<__vector
float>(vec_splat_u32(0)) };
141 static inline SimdFInt32 gmx_simdcall
simdLoad(const std::int32_t* m
, SimdFInt32Tag
)
143 return { vec_ld(0, const_cast<int*>(m
)) };
146 static inline void gmx_simdcall
store(std::int32_t* m
, SimdFInt32 a
)
148 vec_st(a
.simdInternal_
, 0, const_cast<int*>(m
));
151 static inline SimdFInt32 gmx_simdcall
setZeroFI()
153 return { vec_splat_s32(0) };
156 static inline SimdFloat gmx_simdcall
operator&(SimdFloat a
, SimdFloat b
)
158 return { vec_and(a
.simdInternal_
, b
.simdInternal_
) };
161 static inline SimdFloat gmx_simdcall
andNot(SimdFloat a
, SimdFloat b
)
163 return { vec_andc(b
.simdInternal_
, a
.simdInternal_
) };
166 static inline SimdFloat gmx_simdcall
operator|(SimdFloat a
, SimdFloat b
)
168 return { vec_or(a
.simdInternal_
, b
.simdInternal_
) };
171 static inline SimdFloat gmx_simdcall
operator^(SimdFloat a
, SimdFloat b
)
173 return { vec_xor(a
.simdInternal_
, b
.simdInternal_
) };
176 static inline SimdFloat gmx_simdcall
operator+(SimdFloat a
, SimdFloat b
)
178 return { vec_add(a
.simdInternal_
, b
.simdInternal_
) };
181 static inline SimdFloat gmx_simdcall
operator-(SimdFloat a
, SimdFloat b
)
183 return { vec_sub(a
.simdInternal_
, b
.simdInternal_
) };
186 static inline SimdFloat gmx_simdcall
operator-(SimdFloat x
)
188 return { vec_xor(x
.simdInternal_
,
189 reinterpret_cast<__vector
float>(vec_sl(vec_splat_u32(-1), vec_splat_u32(-1)))) };
192 static inline SimdFloat gmx_simdcall
operator*(SimdFloat a
, SimdFloat b
)
194 return { vec_madd(a
.simdInternal_
, b
.simdInternal_
,
195 reinterpret_cast<__vector
float>(vec_splat_u32(0))) };
198 static inline SimdFloat gmx_simdcall
fma(SimdFloat a
, SimdFloat b
, SimdFloat c
)
200 return { vec_madd(a
.simdInternal_
, b
.simdInternal_
, c
.simdInternal_
) };
203 static inline SimdFloat gmx_simdcall
fms(SimdFloat a
, SimdFloat b
, SimdFloat c
)
205 return { vec_madd(a
.simdInternal_
, b
.simdInternal_
, -c
.simdInternal_
) };
208 static inline SimdFloat gmx_simdcall
fnma(SimdFloat a
, SimdFloat b
, SimdFloat c
)
210 return { vec_nmsub(a
.simdInternal_
, b
.simdInternal_
, c
.simdInternal_
) };
213 static inline SimdFloat gmx_simdcall
fnms(SimdFloat a
, SimdFloat b
, SimdFloat c
)
215 return { -vec_madd(a
.simdInternal_
, b
.simdInternal_
, c
.simdInternal_
) };
218 static inline SimdFloat gmx_simdcall
rsqrt(SimdFloat x
)
220 return { vec_rsqrte(x
.simdInternal_
) };
223 static inline SimdFloat gmx_simdcall
rcp(SimdFloat x
)
225 return { vec_re(x
.simdInternal_
) };
228 static inline SimdFloat gmx_simdcall
maskAdd(SimdFloat a
, SimdFloat b
, SimdFBool m
)
230 return { vec_add(a
.simdInternal_
,
231 vec_and(b
.simdInternal_
, reinterpret_cast<__vector
float>(m
.simdInternal_
))) };
234 static inline SimdFloat gmx_simdcall
maskzMul(SimdFloat a
, SimdFloat b
, SimdFBool m
)
236 SimdFloat prod
= a
* b
;
238 return { vec_and(prod
.simdInternal_
, reinterpret_cast<__vector
float>(m
.simdInternal_
)) };
241 static inline SimdFloat gmx_simdcall
maskzFma(SimdFloat a
, SimdFloat b
, SimdFloat c
, SimdFBool m
)
243 SimdFloat prod
= fma(a
, b
, c
);
245 return { vec_and(prod
.simdInternal_
, reinterpret_cast<__vector
float>(m
.simdInternal_
)) };
248 static inline SimdFloat gmx_simdcall
maskzRsqrt(SimdFloat x
, SimdFBool m
)
252 x
.simdInternal_
= vec_sel(one
.simdInternal_
, x
.simdInternal_
, m
.simdInternal_
);
254 return { vec_and(vec_rsqrte(x
.simdInternal_
), reinterpret_cast<__vector
float>(m
.simdInternal_
)) };
257 static inline SimdFloat gmx_simdcall
maskzRcp(SimdFloat x
, SimdFBool m
)
261 x
.simdInternal_
= vec_sel(one
.simdInternal_
, x
.simdInternal_
, m
.simdInternal_
);
263 return { vec_and(vec_re(x
.simdInternal_
), reinterpret_cast<__vector
float>(m
.simdInternal_
)) };
266 static inline SimdFloat gmx_simdcall
abs(SimdFloat x
)
268 return { vec_abs(x
.simdInternal_
) };
271 static inline SimdFloat gmx_simdcall
max(SimdFloat a
, SimdFloat b
)
273 return { vec_max(a
.simdInternal_
, b
.simdInternal_
) };
276 static inline SimdFloat gmx_simdcall
min(SimdFloat a
, SimdFloat b
)
278 return { vec_min(a
.simdInternal_
, b
.simdInternal_
) };
281 static inline SimdFloat gmx_simdcall
round(SimdFloat x
)
283 return { vec_round(x
.simdInternal_
) };
286 static inline SimdFloat gmx_simdcall
trunc(SimdFloat x
)
288 return { vec_trunc(x
.simdInternal_
) };
291 template<MathOptimization opt
= MathOptimization::Safe
>
292 static inline SimdFloat gmx_simdcall
frexp(SimdFloat value
, SimdFInt32
* exponent
)
294 // Generate constants without memory operations
295 const __vector
signed int exponentMask
=
296 vec_sl(vec_add(vec_splat_s32(15), vec_sl(vec_splat_s32(15), vec_splat_u32(4))),
297 vec_add(vec_splat_u32(15), vec_splat_u32(8))); // 0x7F800000
298 const __vector
signed int exponentBias
=
299 vec_sub(vec_sl(vec_splat_s32(1), vec_splat_u32(7)), vec_splat_s32(2)); // 126
300 const SimdFloat
half(0.5F
);
301 __vector
signed int iExponent
;
303 __vector vmxBool
int valueIsZero
=
304 vec_cmpeq(value
.simdInternal_
, reinterpret_cast<__vector
float>(vec_splat_u32(0)));
306 iExponent
= vec_and(reinterpret_cast<__vector
signed int>(value
.simdInternal_
), exponentMask
);
307 iExponent
= vec_sr(iExponent
, vec_add(vec_splat_u32(15), vec_splat_u32(8)));
308 iExponent
= vec_sub(iExponent
, exponentBias
);
309 iExponent
= vec_andc(iExponent
, reinterpret_cast<__vector
int>(valueIsZero
));
311 __vector
float result
=
312 vec_or(vec_andc(value
.simdInternal_
, reinterpret_cast<__vector
float>(exponentMask
)),
314 result
= vec_sel(result
, value
.simdInternal_
, valueIsZero
);
316 exponent
->simdInternal_
= iExponent
;
321 template<MathOptimization opt
= MathOptimization::Safe
>
322 static inline SimdFloat gmx_simdcall
ldexp(SimdFloat value
, SimdFInt32 exponent
)
324 const __vector
signed int exponentBias
=
325 vec_sub(vec_sl(vec_splat_s32(1), vec_splat_u32(7)), vec_splat_s32(1)); // 127
326 __vector
signed int iExponent
;
328 iExponent
= vec_add(exponent
.simdInternal_
, exponentBias
);
330 if (opt
== MathOptimization::Safe
)
332 // Make sure biased argument is not negative
333 iExponent
= vec_max(iExponent
, vec_splat_s32(0));
336 iExponent
= vec_sl(iExponent
, vec_add(vec_splat_u32(15), vec_splat_u32(8)));
338 return { vec_madd(value
.simdInternal_
, reinterpret_cast<__vector
float>(iExponent
),
339 reinterpret_cast<__vector
float>(vec_splat_u32(0))) };
342 static inline float gmx_simdcall
reduce(SimdFloat a
)
344 __vector
float c
= a
.simdInternal_
;
348 c
= vec_add(c
, vec_sld(c
, c
, 8));
349 c
= vec_add(c
, vec_sld(c
, c
, 4));
354 static inline SimdFBool gmx_simdcall
operator==(SimdFloat a
, SimdFloat b
)
356 return { vec_cmpeq(a
.simdInternal_
, b
.simdInternal_
) };
359 static inline SimdFBool gmx_simdcall
operator!=(SimdFloat a
, SimdFloat b
)
361 return { vec_or(vec_cmpgt(a
.simdInternal_
, b
.simdInternal_
),
362 vec_cmplt(a
.simdInternal_
, b
.simdInternal_
)) };
365 static inline SimdFBool gmx_simdcall
operator<(SimdFloat a
, SimdFloat b
)
367 return { vec_cmplt(a
.simdInternal_
, b
.simdInternal_
) };
370 static inline SimdFBool gmx_simdcall
operator<=(SimdFloat a
, SimdFloat b
)
372 return { vec_cmple(a
.simdInternal_
, b
.simdInternal_
) };
375 static inline SimdFBool gmx_simdcall
testBits(SimdFloat a
)
377 return { vec_cmpgt(reinterpret_cast<__vector
unsigned int>(a
.simdInternal_
), vec_splat_u32(0)) };
380 static inline SimdFBool gmx_simdcall
operator&&(SimdFBool a
, SimdFBool b
)
382 return { vec_and(a
.simdInternal_
, b
.simdInternal_
) };
385 static inline SimdFBool gmx_simdcall
operator||(SimdFBool a
, SimdFBool b
)
387 return { vec_or(a
.simdInternal_
, b
.simdInternal_
) };
390 static inline bool gmx_simdcall
anyTrue(SimdFBool a
)
392 return vec_any_ne(a
.simdInternal_
, reinterpret_cast<__vector vmxBool
int>(vec_splat_u32(0)));
395 static inline SimdFloat gmx_simdcall
selectByMask(SimdFloat a
, SimdFBool m
)
397 return { vec_and(a
.simdInternal_
, reinterpret_cast<__vector
float>(m
.simdInternal_
)) };
400 static inline SimdFloat gmx_simdcall
selectByNotMask(SimdFloat a
, SimdFBool m
)
402 return { vec_andc(a
.simdInternal_
, reinterpret_cast<__vector
float>(m
.simdInternal_
)) };
405 static inline SimdFloat gmx_simdcall
blend(SimdFloat a
, SimdFloat b
, SimdFBool sel
)
407 return { vec_sel(a
.simdInternal_
, b
.simdInternal_
, sel
.simdInternal_
) };
410 static inline SimdFInt32 gmx_simdcall
operator&(SimdFInt32 a
, SimdFInt32 b
)
412 return { vec_and(a
.simdInternal_
, b
.simdInternal_
) };
415 static inline SimdFInt32 gmx_simdcall
andNot(SimdFInt32 a
, SimdFInt32 b
)
417 return { vec_andc(b
.simdInternal_
, a
.simdInternal_
) };
420 static inline SimdFInt32 gmx_simdcall
operator|(SimdFInt32 a
, SimdFInt32 b
)
422 return { vec_or(a
.simdInternal_
, b
.simdInternal_
) };
425 static inline SimdFInt32 gmx_simdcall
operator^(SimdFInt32 a
, SimdFInt32 b
)
427 return { vec_xor(a
.simdInternal_
, b
.simdInternal_
) };
430 static inline SimdFInt32 gmx_simdcall
operator+(SimdFInt32 a
, SimdFInt32 b
)
432 return { vec_add(a
.simdInternal_
, b
.simdInternal_
) };
435 static inline SimdFInt32 gmx_simdcall
operator-(SimdFInt32 a
, SimdFInt32 b
)
437 return { vec_sub(a
.simdInternal_
, b
.simdInternal_
) };
440 static inline SimdFInt32 gmx_simdcall
operator*(SimdFInt32 a
, SimdFInt32 b
)
442 return { a
.simdInternal_
* b
.simdInternal_
};
445 static inline SimdFIBool gmx_simdcall
operator==(SimdFInt32 a
, SimdFInt32 b
)
447 return { vec_cmpeq(a
.simdInternal_
, b
.simdInternal_
) };
450 static inline SimdFIBool gmx_simdcall
testBits(SimdFInt32 a
)
452 return { vec_cmpgt(reinterpret_cast<__vector
unsigned int>(a
.simdInternal_
), vec_splat_u32(0)) };
455 static inline SimdFIBool gmx_simdcall
operator<(SimdFInt32 a
, SimdFInt32 b
)
457 return { vec_cmplt(a
.simdInternal_
, b
.simdInternal_
) };
460 static inline SimdFIBool gmx_simdcall
operator&&(SimdFIBool a
, SimdFIBool b
)
462 return { vec_and(a
.simdInternal_
, b
.simdInternal_
) };
465 static inline SimdFIBool gmx_simdcall
operator||(SimdFIBool a
, SimdFIBool b
)
467 return { vec_or(a
.simdInternal_
, b
.simdInternal_
) };
470 static inline bool gmx_simdcall
anyTrue(SimdFIBool a
)
472 return vec_any_ne(a
.simdInternal_
, reinterpret_cast<__vector vmxBool
int>(vec_splat_u32(0)));
475 static inline SimdFInt32 gmx_simdcall
selectByMask(SimdFInt32 a
, SimdFIBool m
)
477 return { vec_and(a
.simdInternal_
, reinterpret_cast<__vector
signed int>(m
.simdInternal_
)) };
480 static inline SimdFInt32 gmx_simdcall
selectByNotMask(SimdFInt32 a
, SimdFIBool m
)
482 return { vec_andc(a
.simdInternal_
, reinterpret_cast<__vector
signed int>(m
.simdInternal_
)) };
485 static inline SimdFInt32 gmx_simdcall
blend(SimdFInt32 a
, SimdFInt32 b
, SimdFIBool sel
)
487 return { vec_sel(a
.simdInternal_
, b
.simdInternal_
, sel
.simdInternal_
) };
490 static inline SimdFInt32 gmx_simdcall
cvtR2I(SimdFloat a
)
492 return { vec_cts(vec_round(a
.simdInternal_
), 0) };
495 static inline SimdFInt32 gmx_simdcall
cvttR2I(SimdFloat a
)
497 return { vec_cts(a
.simdInternal_
, 0) };
500 static inline SimdFloat gmx_simdcall
cvtI2R(SimdFInt32 a
)
502 return { vec_ctf(a
.simdInternal_
, 0) };
505 static inline SimdFIBool gmx_simdcall
cvtB2IB(SimdFBool a
)
507 return { a
.simdInternal_
};
510 static inline SimdFBool gmx_simdcall
cvtIB2B(SimdFIBool a
)
512 return { a
.simdInternal_
};
517 #endif // GMX_SIMD_IMPLEMENTATION_IBM_VMX_SIMD_FLOAT_H