2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2014,2015,2017, 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_QPX_SIMD4_DOUBLE_H
37 #define GMX_SIMD_IMPLEMENTATION_IBM_QPX_SIMD4_DOUBLE_H
43 // Assert is buggy on xlc with high optimization, so we skip it for QPX
58 Simd4Double(double d) : simdInternal_(vec_splats(d)) {}
60 // Internal utility constructor to simplify return statements
61 Simd4Double(vector4double simd) : simdInternal_(simd) {}
63 vector4double simdInternal_;
71 //! \brief Construct from scalar bool
72 Simd4DBool(bool b) : simdInternal_(vec_splats(b ? 1.0 : -1.0)) {}
74 // Internal utility constructor to simplify return statements
75 Simd4DBool(vector4double simd) : simdInternal_(simd) {}
77 vector4double simdInternal_;
80 static inline Simd4Double gmx_simdcall
81 load4(const double *m)
85 vec_ld(0, const_cast<double *>(m))
89 vec_lda(0, const_cast<double *>(m))
94 static inline void gmx_simdcall
95 store4(double *m, Simd4Double a)
98 vec_st(a.simdInternal_, 0, m);
100 vec_sta(a.simdInternal_, 0, m);
104 static inline Simd4Double gmx_simdcall
112 static inline Simd4Double gmx_simdcall
113 operator+(Simd4Double a, Simd4Double b)
116 vec_add(a.simdInternal_, b.simdInternal_)
120 static inline Simd4Double gmx_simdcall
121 operator-(Simd4Double a, Simd4Double b)
124 vec_sub(a.simdInternal_, b.simdInternal_)
128 static inline Simd4Double gmx_simdcall
129 operator-(Simd4Double x)
132 vec_neg(x.simdInternal_)
136 static inline Simd4Double gmx_simdcall
137 operator*(Simd4Double a, Simd4Double b)
140 vec_mul(a.simdInternal_, b.simdInternal_)
144 static inline Simd4Double gmx_simdcall
145 fma(Simd4Double a, Simd4Double b, Simd4Double c)
148 vec_madd(a.simdInternal_, b.simdInternal_, c.simdInternal_)
152 static inline Simd4Double gmx_simdcall
153 fms(Simd4Double a, Simd4Double b, Simd4Double c)
156 vec_msub(a.simdInternal_, b.simdInternal_, c.simdInternal_)
160 static inline Simd4Double gmx_simdcall
161 fnma(Simd4Double a, Simd4Double b, Simd4Double c)
164 vec_nmsub(a.simdInternal_, b.simdInternal_, c.simdInternal_)
168 static inline Simd4Double gmx_simdcall
169 fnms(Simd4Double a, Simd4Double b, Simd4Double c)
172 vec_nmadd(a.simdInternal_, b.simdInternal_, c.simdInternal_)
176 static inline Simd4Double gmx_simdcall
180 vec_rsqrte(x.simdInternal_)
184 static inline Simd4Double gmx_simdcall
188 vec_abs( x.simdInternal_ )
192 static inline Simd4Double gmx_simdcall
193 max(Simd4Double a, Simd4Double b)
196 vec_sel(b.simdInternal_, a.simdInternal_, vec_sub(a.simdInternal_, b.simdInternal_))
200 static inline Simd4Double gmx_simdcall
201 min(Simd4Double a, Simd4Double b)
204 vec_sel(b.simdInternal_, a.simdInternal_, vec_sub(b.simdInternal_, a.simdInternal_))
208 static inline Simd4Double gmx_simdcall
211 // Note: It is critical to use vec_cfid(vec_ctid(a)) for the implementation
212 // here, since vec_round() does not adhere to the FP control
213 // word rounding scheme. We rely on float-to-float and float-to-integer
214 // rounding being the same for half-way values in a few algorithms.
217 vec_cfid(vec_ctid(x.simdInternal_))
221 static inline Simd4Double gmx_simdcall
225 vec_trunc(x.simdInternal_)
229 static inline double gmx_simdcall
230 dotProduct(Simd4Double a, Simd4Double b)
232 vector4double dp_sh0 = vec_mul(a.simdInternal_, b.simdInternal_);
233 vector4double dp_sh1 = vec_sldw(dp_sh0, dp_sh0, 1);
234 vector4double dp_sh2 = vec_sldw(dp_sh0, dp_sh0, 2);
235 vector4double dp = vec_add(dp_sh2, vec_add(dp_sh0, dp_sh1));
237 return vec_extract(dp, 0);
240 static inline void gmx_simdcall
241 transpose(Simd4Double * v0, Simd4Double * v1,
242 Simd4Double * v2, Simd4Double * v3)
244 vector4double t0 = vec_perm(v0->simdInternal_, v2->simdInternal_, vec_gpci(00415));
245 vector4double t1 = vec_perm(v0->simdInternal_, v2->simdInternal_, vec_gpci(02637));
246 vector4double t2 = vec_perm(v1->simdInternal_, v3->simdInternal_, vec_gpci(00415));
247 vector4double t3 = vec_perm(v1->simdInternal_, v3->simdInternal_, vec_gpci(02637));
248 v0->simdInternal_ = vec_perm(t0, t2, vec_gpci(00415));
249 v1->simdInternal_ = vec_perm(t0, t2, vec_gpci(02637));
250 v2->simdInternal_ = vec_perm(t1, t3, vec_gpci(00415));
251 v3->simdInternal_ = vec_perm(t1, t3, vec_gpci(02637));
254 static inline Simd4DBool gmx_simdcall
255 operator==(Simd4Double a, Simd4Double b)
258 vec_cmpeq(a.simdInternal_, b.simdInternal_)
262 static inline Simd4DBool gmx_simdcall
263 operator!=(Simd4Double a, Simd4Double b)
266 vec_not(vec_cmpeq(a.simdInternal_, b.simdInternal_))
270 static inline Simd4DBool gmx_simdcall
271 operator<(Simd4Double a, Simd4Double b)
274 vec_cmplt(a.simdInternal_, b.simdInternal_)
278 static inline Simd4DBool gmx_simdcall
279 operator<=(Simd4Double a, Simd4Double b)
282 vec_or(vec_cmplt(a.simdInternal_, b.simdInternal_), vec_cmpeq(a.simdInternal_, b.simdInternal_))
286 static inline Simd4DBool gmx_simdcall
287 operator&&(Simd4DBool a, Simd4DBool b)
290 vec_and(a.simdInternal_, b.simdInternal_)
294 static inline Simd4DBool gmx_simdcall
295 operator||(Simd4DBool a, Simd4DBool b)
298 vec_or(a.simdInternal_, b.simdInternal_)
302 static inline bool gmx_simdcall
303 anyTrue(Simd4DBool a)
305 vector4double b = vec_sldw(a.simdInternal_, a.simdInternal_, 2);
307 a.simdInternal_ = vec_or(a.simdInternal_, b);
308 b = vec_sldw(a.simdInternal_, a.simdInternal_, 1);
309 b = vec_or(a.simdInternal_, b);
310 return (vec_extract(b, 0) > 0);
313 static inline Simd4Double gmx_simdcall
314 selectByMask(Simd4Double a, Simd4DBool m)
317 vec_sel(vec_splats(0.0), a.simdInternal_, m.simdInternal_)
321 static inline Simd4Double gmx_simdcall
322 selectByNotMask(Simd4Double a, Simd4DBool m)
325 vec_sel(a.simdInternal_, vec_splats(0.0), m.simdInternal_)
329 static inline Simd4Double gmx_simdcall
330 blend(Simd4Double a, Simd4Double b, Simd4DBool sel)
333 vec_sel(a.simdInternal_, b.simdInternal_, sel.simdInternal_)
337 static inline double gmx_simdcall
338 reduce(Simd4Double x)
340 vector4double y = vec_sldw(x.simdInternal_, x.simdInternal_, 2);
343 y = vec_add(y, x.simdInternal_);
344 z = vec_sldw(y, y, 1);
346 return vec_extract(y, 0);
351 #endif // GMX_SIMD_IMPLEMENTATION_IBM_QPX_SIMD4_DOUBLE_H