Apply re-formatting to C++ in src/ tree.
[alexxy/gromacs.git] / src / gromacs / simd / impl_x86_avx_256 / impl_x86_avx_256_simd_double.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2014,2015,2016,2017,2018 by the GROMACS development team.
5  * Copyright (c) 2019,2020, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36
37 #ifndef GMX_SIMD_IMPL_X86_AVX_256_SIMD_DOUBLE_H
38 #define GMX_SIMD_IMPL_X86_AVX_256_SIMD_DOUBLE_H
39
40 #include "config.h"
41
42 #include <cassert>
43 #include <cstddef>
44 #include <cstdint>
45
46 #include <immintrin.h>
47
48 #include "gromacs/math/utilities.h"
49
50 #include "impl_x86_avx_256_simd_float.h"
51
52
53 namespace gmx
54 {
55
56 class SimdDouble
57 {
58 public:
59     MSVC_DIAGNOSTIC_IGNORE(26495) // simdInternal_ is not being initialized!
60     SimdDouble() {}
61     MSVC_DIAGNOSTIC_RESET
62     SimdDouble(double d) : simdInternal_(_mm256_set1_pd(d)) {}
63
64     // Internal utility constructor to simplify return statements
65     SimdDouble(__m256d simd) : simdInternal_(simd) {}
66
67     __m256d simdInternal_;
68 };
69
70 class SimdDInt32
71 {
72 public:
73     MSVC_DIAGNOSTIC_IGNORE(26495) // simdInternal_ is not being initialized!
74     SimdDInt32() {}
75     MSVC_DIAGNOSTIC_RESET
76     SimdDInt32(std::int32_t i) : simdInternal_(_mm_set1_epi32(i)) {}
77
78     // Internal utility constructor to simplify return statements
79     SimdDInt32(__m128i simd) : simdInternal_(simd) {}
80
81     __m128i simdInternal_;
82 };
83
84 class SimdDBool
85 {
86 public:
87     MSVC_DIAGNOSTIC_IGNORE(26495) // simdInternal_ is not being initialized!
88     SimdDBool() {}
89     MSVC_DIAGNOSTIC_RESET
90     SimdDBool(bool b) : simdInternal_(_mm256_castsi256_pd(_mm256_set1_epi32(b ? 0xFFFFFFFF : 0))) {}
91
92     // Internal utility constructor to simplify return statements
93     SimdDBool(__m256d simd) : simdInternal_(simd) {}
94
95     __m256d simdInternal_;
96 };
97
98 class SimdDIBool
99 {
100 public:
101     MSVC_DIAGNOSTIC_IGNORE(26495) // simdInternal_ is not being initialized!
102     SimdDIBool() {}
103     MSVC_DIAGNOSTIC_RESET
104     SimdDIBool(bool b) : simdInternal_(_mm_set1_epi32(b ? 0xFFFFFFFF : 0)) {}
105
106     // Internal utility constructor to simplify return statements
107     SimdDIBool(__m128i simd) : simdInternal_(simd) {}
108
109     __m128i simdInternal_;
110 };
111
112
113 static inline SimdDouble gmx_simdcall simdLoad(const double* m, SimdDoubleTag /*unused*/ = {})
114 {
115     assert(std::size_t(m) % 32 == 0);
116     return { _mm256_load_pd(m) };
117 }
118
119 static inline void gmx_simdcall store(double* m, SimdDouble a)
120 {
121     assert(std::size_t(m) % 32 == 0);
122     _mm256_store_pd(m, a.simdInternal_);
123 }
124
125 static inline SimdDouble gmx_simdcall simdLoadU(const double* m, SimdDoubleTag /*unused*/ = {})
126 {
127     return { _mm256_loadu_pd(m) };
128 }
129
130 static inline void gmx_simdcall storeU(double* m, SimdDouble a)
131 {
132     _mm256_storeu_pd(m, a.simdInternal_);
133 }
134
135 static inline SimdDouble gmx_simdcall setZeroD()
136 {
137     return { _mm256_setzero_pd() };
138 }
139
140 static inline SimdDInt32 gmx_simdcall simdLoad(const std::int32_t* m, SimdDInt32Tag /*unused*/)
141 {
142     assert(std::size_t(m) % 16 == 0);
143     return { _mm_load_si128(reinterpret_cast<const __m128i*>(m)) };
144 }
145
146 static inline void gmx_simdcall store(std::int32_t* m, SimdDInt32 a)
147 {
148     assert(std::size_t(m) % 16 == 0);
149     _mm_store_si128(reinterpret_cast<__m128i*>(m), a.simdInternal_);
150 }
151
152 static inline SimdDInt32 gmx_simdcall simdLoadU(const std::int32_t* m, SimdDInt32Tag /*unused*/)
153 {
154     return { _mm_loadu_si128(reinterpret_cast<const __m128i*>(m)) };
155 }
156
157 static inline void gmx_simdcall storeU(std::int32_t* m, SimdDInt32 a)
158 {
159     _mm_storeu_si128(reinterpret_cast<__m128i*>(m), a.simdInternal_);
160 }
161
162 static inline SimdDInt32 gmx_simdcall setZeroDI()
163 {
164     return { _mm_setzero_si128() };
165 }
166
167 template<int index>
168 static inline std::int32_t gmx_simdcall extract(SimdDInt32 a)
169 {
170     return _mm_extract_epi32(a.simdInternal_, index);
171 }
172
173 static inline SimdDouble gmx_simdcall operator&(SimdDouble a, SimdDouble b)
174 {
175     return { _mm256_and_pd(a.simdInternal_, b.simdInternal_) };
176 }
177
178 static inline SimdDouble gmx_simdcall andNot(SimdDouble a, SimdDouble b)
179 {
180     return { _mm256_andnot_pd(a.simdInternal_, b.simdInternal_) };
181 }
182
183 static inline SimdDouble gmx_simdcall operator|(SimdDouble a, SimdDouble b)
184 {
185     return { _mm256_or_pd(a.simdInternal_, b.simdInternal_) };
186 }
187
188 static inline SimdDouble gmx_simdcall operator^(SimdDouble a, SimdDouble b)
189 {
190     return { _mm256_xor_pd(a.simdInternal_, b.simdInternal_) };
191 }
192
193 static inline SimdDouble gmx_simdcall operator+(SimdDouble a, SimdDouble b)
194 {
195     return { _mm256_add_pd(a.simdInternal_, b.simdInternal_) };
196 }
197
198 static inline SimdDouble gmx_simdcall operator-(SimdDouble a, SimdDouble b)
199 {
200     return { _mm256_sub_pd(a.simdInternal_, b.simdInternal_) };
201 }
202
203 static inline SimdDouble gmx_simdcall operator-(SimdDouble x)
204 {
205     return { _mm256_xor_pd(x.simdInternal_, _mm256_set1_pd(GMX_DOUBLE_NEGZERO)) };
206 }
207
208 static inline SimdDouble gmx_simdcall operator*(SimdDouble a, SimdDouble b)
209 {
210     return { _mm256_mul_pd(a.simdInternal_, b.simdInternal_) };
211 }
212
213 // Override for AVX2 and higher
214 #if GMX_SIMD_X86_AVX_256
215 static inline SimdDouble gmx_simdcall fma(SimdDouble a, SimdDouble b, SimdDouble c)
216 {
217     return { _mm256_add_pd(_mm256_mul_pd(a.simdInternal_, b.simdInternal_), c.simdInternal_) };
218 }
219
220 static inline SimdDouble gmx_simdcall fms(SimdDouble a, SimdDouble b, SimdDouble c)
221 {
222     return { _mm256_sub_pd(_mm256_mul_pd(a.simdInternal_, b.simdInternal_), c.simdInternal_) };
223 }
224
225 static inline SimdDouble gmx_simdcall fnma(SimdDouble a, SimdDouble b, SimdDouble c)
226 {
227     return { _mm256_sub_pd(c.simdInternal_, _mm256_mul_pd(a.simdInternal_, b.simdInternal_)) };
228 }
229
230 static inline SimdDouble gmx_simdcall fnms(SimdDouble a, SimdDouble b, SimdDouble c)
231 {
232     return { _mm256_sub_pd(_mm256_setzero_pd(),
233                            _mm256_add_pd(_mm256_mul_pd(a.simdInternal_, b.simdInternal_), c.simdInternal_)) };
234 }
235 #endif
236
237 static inline SimdDouble gmx_simdcall rsqrt(SimdDouble x)
238 {
239     return { _mm256_cvtps_pd(_mm_rsqrt_ps(_mm256_cvtpd_ps(x.simdInternal_))) };
240 }
241
242 static inline SimdDouble gmx_simdcall rcp(SimdDouble x)
243 {
244     return { _mm256_cvtps_pd(_mm_rcp_ps(_mm256_cvtpd_ps(x.simdInternal_))) };
245 }
246
247 static inline SimdDouble gmx_simdcall maskAdd(SimdDouble a, SimdDouble b, SimdDBool m)
248 {
249     return { _mm256_add_pd(a.simdInternal_, _mm256_and_pd(b.simdInternal_, m.simdInternal_)) };
250 }
251
252 static inline SimdDouble gmx_simdcall maskzMul(SimdDouble a, SimdDouble b, SimdDBool m)
253 {
254     return { _mm256_and_pd(_mm256_mul_pd(a.simdInternal_, b.simdInternal_), m.simdInternal_) };
255 }
256
257 static inline SimdDouble maskzFma(SimdDouble a, SimdDouble b, SimdDouble c, SimdDBool m)
258 {
259     return { _mm256_and_pd(_mm256_add_pd(_mm256_mul_pd(a.simdInternal_, b.simdInternal_), c.simdInternal_),
260                            m.simdInternal_) };
261 }
262
263 static inline SimdDouble maskzRsqrt(SimdDouble x, SimdDBool m)
264 {
265 #ifndef NDEBUG
266     x.simdInternal_ = _mm256_blendv_pd(_mm256_set1_pd(1.0), x.simdInternal_, m.simdInternal_);
267 #endif
268     return { _mm256_and_pd(_mm256_cvtps_pd(_mm_rsqrt_ps(_mm256_cvtpd_ps(x.simdInternal_))), m.simdInternal_) };
269 }
270
271 static inline SimdDouble maskzRcp(SimdDouble x, SimdDBool m)
272 {
273 #ifndef NDEBUG
274     x.simdInternal_ = _mm256_blendv_pd(_mm256_set1_pd(1.0), x.simdInternal_, m.simdInternal_);
275 #endif
276     return { _mm256_and_pd(_mm256_cvtps_pd(_mm_rcp_ps(_mm256_cvtpd_ps(x.simdInternal_))), m.simdInternal_) };
277 }
278
279 static inline SimdDouble gmx_simdcall abs(SimdDouble x)
280 {
281     return { _mm256_andnot_pd(_mm256_set1_pd(GMX_DOUBLE_NEGZERO), x.simdInternal_) };
282 }
283
284 static inline SimdDouble gmx_simdcall max(SimdDouble a, SimdDouble b)
285 {
286     return { _mm256_max_pd(a.simdInternal_, b.simdInternal_) };
287 }
288
289 static inline SimdDouble gmx_simdcall min(SimdDouble a, SimdDouble b)
290 {
291     return { _mm256_min_pd(a.simdInternal_, b.simdInternal_) };
292 }
293
294 static inline SimdDouble gmx_simdcall round(SimdDouble x)
295 {
296     return { _mm256_round_pd(x.simdInternal_, _MM_FROUND_NINT) };
297 }
298
299 static inline SimdDouble gmx_simdcall trunc(SimdDouble x)
300 {
301     return { _mm256_round_pd(x.simdInternal_, _MM_FROUND_TRUNC) };
302 }
303
304 // Override for AVX2 and higher
305 #if GMX_SIMD_X86_AVX_256
306 template<MathOptimization opt = MathOptimization::Safe>
307 static inline SimdDouble frexp(SimdDouble value, SimdDInt32* exponent)
308 {
309     const __m256d exponentMask = _mm256_castsi256_pd(_mm256_set1_epi64x(0x7FF0000000000000LL));
310     const __m256d mantissaMask = _mm256_castsi256_pd(_mm256_set1_epi64x(0x800FFFFFFFFFFFFFLL));
311     const __m256d half         = _mm256_set1_pd(0.5);
312     const __m128i exponentBias = _mm_set1_epi32(1022); // add 1 to make our definition identical to frexp()
313
314     __m256i iExponent     = _mm256_castpd_si256(_mm256_and_pd(value.simdInternal_, exponentMask));
315     __m128i iExponentHigh = _mm256_extractf128_si256(iExponent, 0x1);
316     __m128i iExponentLow  = _mm256_castsi256_si128(iExponent);
317     iExponentLow          = _mm_srli_epi64(iExponentLow, 52);
318     iExponentHigh         = _mm_srli_epi64(iExponentHigh, 52);
319     iExponentLow          = _mm_shuffle_epi32(iExponentLow, _MM_SHUFFLE(1, 1, 2, 0));
320     iExponentHigh         = _mm_shuffle_epi32(iExponentHigh, _MM_SHUFFLE(2, 0, 1, 1));
321     // We need to store the return in a 128-bit integer variable, so reuse iExponentLow for both
322     iExponentLow = _mm_or_si128(iExponentLow, iExponentHigh);
323     iExponentLow = _mm_sub_epi32(iExponentLow, exponentBias);
324
325     __m256d result = _mm256_or_pd(_mm256_and_pd(value.simdInternal_, mantissaMask), half);
326
327     if (opt == MathOptimization::Safe)
328     {
329         __m256d valueIsZero = _mm256_cmp_pd(_mm256_setzero_pd(), value.simdInternal_, _CMP_EQ_OQ);
330         // This looks more complex than it is: the valueIsZero variable contains 4x 64-bit double
331         // precision fields, but a bit below we'll need a corresponding integer variable with 4x
332         // 32-bit fields. Since AVX1 does not support shuffling across the upper/lower 128-bit
333         // lanes, we need to extract those first, and then shuffle between two 128-bit variables.
334         __m128i iValueIsZero = _mm_castps_si128(
335                 _mm_shuffle_ps(_mm256_extractf128_ps(_mm256_castpd_ps(valueIsZero), 0x0),
336                                _mm256_extractf128_ps(_mm256_castpd_ps(valueIsZero), 0x1),
337                                _MM_SHUFFLE(2, 0, 2, 0)));
338
339         // Set exponent to 0 when input value was zero
340         iExponentLow = _mm_andnot_si128(iValueIsZero, iExponentLow);
341
342         // Set result to +-0 if the corresponding input value was +-0
343         result = _mm256_blendv_pd(result, value.simdInternal_, valueIsZero);
344     }
345     exponent->simdInternal_ = iExponentLow;
346
347     return { result };
348 }
349
350 template<MathOptimization opt = MathOptimization::Safe>
351 static inline SimdDouble ldexp(SimdDouble value, SimdDInt32 exponent)
352 {
353     const __m128i exponentBias = _mm_set1_epi32(1023);
354     __m128i       iExponentLow, iExponentHigh;
355     __m256d       fExponent;
356
357     iExponentLow = _mm_add_epi32(exponent.simdInternal_, exponentBias);
358
359     if (opt == MathOptimization::Safe)
360     {
361         // Make sure biased argument is not negative
362         iExponentLow = _mm_max_epi32(iExponentLow, _mm_setzero_si128());
363     }
364
365     iExponentHigh = _mm_shuffle_epi32(iExponentLow, _MM_SHUFFLE(3, 3, 2, 2));
366     iExponentLow  = _mm_shuffle_epi32(iExponentLow, _MM_SHUFFLE(1, 1, 0, 0));
367     iExponentHigh = _mm_slli_epi64(iExponentHigh, 52);
368     iExponentLow  = _mm_slli_epi64(iExponentLow, 52);
369     fExponent     = _mm256_castsi256_pd(
370             _mm256_insertf128_si256(_mm256_castsi128_si256(iExponentLow), iExponentHigh, 0x1));
371     return { _mm256_mul_pd(value.simdInternal_, fExponent) };
372 }
373 #endif
374
375 static inline double gmx_simdcall reduce(SimdDouble a)
376 {
377     __m128d a0, a1;
378     a.simdInternal_ = _mm256_add_pd(a.simdInternal_, _mm256_permute_pd(a.simdInternal_, 0b0101));
379     a0              = _mm256_castpd256_pd128(a.simdInternal_);
380     a1              = _mm256_extractf128_pd(a.simdInternal_, 0x1);
381     a0              = _mm_add_sd(a0, a1);
382
383     return *reinterpret_cast<double*>(&a0);
384 }
385
386 static inline SimdDBool gmx_simdcall operator==(SimdDouble a, SimdDouble b)
387 {
388     return { _mm256_cmp_pd(a.simdInternal_, b.simdInternal_, _CMP_EQ_OQ) };
389 }
390
391 static inline SimdDBool gmx_simdcall operator!=(SimdDouble a, SimdDouble b)
392 {
393     return { _mm256_cmp_pd(a.simdInternal_, b.simdInternal_, _CMP_NEQ_OQ) };
394 }
395
396 static inline SimdDBool gmx_simdcall operator<(SimdDouble a, SimdDouble b)
397 {
398     return { _mm256_cmp_pd(a.simdInternal_, b.simdInternal_, _CMP_LT_OQ) };
399 }
400
401 static inline SimdDBool gmx_simdcall operator<=(SimdDouble a, SimdDouble b)
402 {
403     return { _mm256_cmp_pd(a.simdInternal_, b.simdInternal_, _CMP_LE_OQ) };
404 }
405
406 // Override for AVX2 and higher
407 #if GMX_SIMD_X86_AVX_256
408 static inline SimdDBool gmx_simdcall testBits(SimdDouble a)
409 {
410     // Do an or of the low/high 32 bits of each double (so the data is replicated),
411     // and then use the same algorithm as we use for single precision.
412     __m256 tst = _mm256_castpd_ps(a.simdInternal_);
413
414     tst = _mm256_or_ps(tst, _mm256_permute_ps(tst, _MM_SHUFFLE(2, 3, 0, 1)));
415     tst = _mm256_cvtepi32_ps(_mm256_castps_si256(tst));
416
417     return { _mm256_castps_pd(_mm256_cmp_ps(tst, _mm256_setzero_ps(), _CMP_NEQ_OQ)) };
418 }
419 #endif
420
421 static inline SimdDBool gmx_simdcall operator&&(SimdDBool a, SimdDBool b)
422 {
423     return { _mm256_and_pd(a.simdInternal_, b.simdInternal_) };
424 }
425
426 static inline SimdDBool gmx_simdcall operator||(SimdDBool a, SimdDBool b)
427 {
428     return { _mm256_or_pd(a.simdInternal_, b.simdInternal_) };
429 }
430
431 static inline bool gmx_simdcall anyTrue(SimdDBool a)
432 {
433     return _mm256_movemask_pd(a.simdInternal_) != 0;
434 }
435
436 static inline SimdDouble gmx_simdcall selectByMask(SimdDouble a, SimdDBool mask)
437 {
438     return { _mm256_and_pd(a.simdInternal_, mask.simdInternal_) };
439 }
440
441 static inline SimdDouble gmx_simdcall selectByNotMask(SimdDouble a, SimdDBool mask)
442 {
443     return { _mm256_andnot_pd(mask.simdInternal_, a.simdInternal_) };
444 }
445
446 static inline SimdDouble gmx_simdcall blend(SimdDouble a, SimdDouble b, SimdDBool sel)
447 {
448     return { _mm256_blendv_pd(a.simdInternal_, b.simdInternal_, sel.simdInternal_) };
449 }
450
451 static inline SimdDInt32 gmx_simdcall operator&(SimdDInt32 a, SimdDInt32 b)
452 {
453     return { _mm_and_si128(a.simdInternal_, b.simdInternal_) };
454 }
455
456 static inline SimdDInt32 gmx_simdcall andNot(SimdDInt32 a, SimdDInt32 b)
457 {
458     return { _mm_andnot_si128(a.simdInternal_, b.simdInternal_) };
459 }
460
461 static inline SimdDInt32 gmx_simdcall operator|(SimdDInt32 a, SimdDInt32 b)
462 {
463     return { _mm_or_si128(a.simdInternal_, b.simdInternal_) };
464 }
465
466 static inline SimdDInt32 gmx_simdcall operator^(SimdDInt32 a, SimdDInt32 b)
467 {
468     return { _mm_xor_si128(a.simdInternal_, b.simdInternal_) };
469 }
470
471 static inline SimdDInt32 gmx_simdcall operator+(SimdDInt32 a, SimdDInt32 b)
472 {
473     return { _mm_add_epi32(a.simdInternal_, b.simdInternal_) };
474 }
475
476 static inline SimdDInt32 gmx_simdcall operator-(SimdDInt32 a, SimdDInt32 b)
477 {
478     return { _mm_sub_epi32(a.simdInternal_, b.simdInternal_) };
479 }
480
481 static inline SimdDInt32 gmx_simdcall operator*(SimdDInt32 a, SimdDInt32 b)
482 {
483     return { _mm_mullo_epi32(a.simdInternal_, b.simdInternal_) };
484 }
485
486 static inline SimdDIBool gmx_simdcall operator==(SimdDInt32 a, SimdDInt32 b)
487 {
488     return { _mm_cmpeq_epi32(a.simdInternal_, b.simdInternal_) };
489 }
490
491 static inline SimdDIBool gmx_simdcall operator<(SimdDInt32 a, SimdDInt32 b)
492 {
493     return { _mm_cmplt_epi32(a.simdInternal_, b.simdInternal_) };
494 }
495
496 static inline SimdDIBool gmx_simdcall testBits(SimdDInt32 a)
497 {
498     __m128i x   = a.simdInternal_;
499     __m128i res = _mm_andnot_si128(_mm_cmpeq_epi32(x, _mm_setzero_si128()), _mm_cmpeq_epi32(x, x));
500
501     return { res };
502 }
503
504 static inline SimdDIBool gmx_simdcall operator&&(SimdDIBool a, SimdDIBool b)
505 {
506     return { _mm_and_si128(a.simdInternal_, b.simdInternal_) };
507 }
508
509 static inline SimdDIBool gmx_simdcall operator||(SimdDIBool a, SimdDIBool b)
510 {
511     return { _mm_or_si128(a.simdInternal_, b.simdInternal_) };
512 }
513
514 static inline bool gmx_simdcall anyTrue(SimdDIBool a)
515 {
516     return _mm_movemask_epi8(a.simdInternal_) != 0;
517 }
518
519 static inline SimdDInt32 gmx_simdcall selectByMask(SimdDInt32 a, SimdDIBool mask)
520 {
521     return { _mm_and_si128(a.simdInternal_, mask.simdInternal_) };
522 }
523
524 static inline SimdDInt32 gmx_simdcall selectByNotMask(SimdDInt32 a, SimdDIBool mask)
525 {
526     return { _mm_andnot_si128(mask.simdInternal_, a.simdInternal_) };
527 }
528
529 static inline SimdDInt32 gmx_simdcall blend(SimdDInt32 a, SimdDInt32 b, SimdDIBool sel)
530 {
531     return { _mm_blendv_epi8(a.simdInternal_, b.simdInternal_, sel.simdInternal_) };
532 }
533
534 static inline SimdDInt32 gmx_simdcall cvtR2I(SimdDouble a)
535 {
536     return { _mm256_cvtpd_epi32(a.simdInternal_) };
537 }
538
539 static inline SimdDInt32 gmx_simdcall cvttR2I(SimdDouble a)
540 {
541     return { _mm256_cvttpd_epi32(a.simdInternal_) };
542 }
543
544 static inline SimdDouble gmx_simdcall cvtI2R(SimdDInt32 a)
545 {
546     return { _mm256_cvtepi32_pd(a.simdInternal_) };
547 }
548
549 static inline SimdDIBool gmx_simdcall cvtB2IB(SimdDBool a)
550 {
551     __m128i a1 = _mm256_extractf128_si256(_mm256_castpd_si256(a.simdInternal_), 0x1);
552     __m128i a0 = _mm256_castsi256_si128(_mm256_castpd_si256(a.simdInternal_));
553     a0         = _mm_shuffle_epi32(a0, _MM_SHUFFLE(2, 0, 2, 0));
554     a1         = _mm_shuffle_epi32(a1, _MM_SHUFFLE(2, 0, 2, 0));
555
556     return { _mm_blend_epi16(a0, a1, 0xF0) };
557 }
558
559 static inline SimdDBool gmx_simdcall cvtIB2B(SimdDIBool a)
560 {
561     __m128d lo = _mm_castsi128_pd(_mm_unpacklo_epi32(a.simdInternal_, a.simdInternal_));
562     __m128d hi = _mm_castsi128_pd(_mm_unpackhi_epi32(a.simdInternal_, a.simdInternal_));
563
564     return { _mm256_insertf128_pd(_mm256_castpd128_pd256(lo), hi, 0x1) };
565 }
566
567 static inline void gmx_simdcall cvtF2DD(SimdFloat f, SimdDouble* d0, SimdDouble* d1)
568 {
569     d0->simdInternal_ = _mm256_cvtps_pd(_mm256_castps256_ps128(f.simdInternal_));
570     d1->simdInternal_ = _mm256_cvtps_pd(_mm256_extractf128_ps(f.simdInternal_, 0x1));
571 }
572
573 static inline SimdFloat gmx_simdcall cvtDD2F(SimdDouble d0, SimdDouble d1)
574 {
575     __m128 f0 = _mm256_cvtpd_ps(d0.simdInternal_);
576     __m128 f1 = _mm256_cvtpd_ps(d1.simdInternal_);
577     return { _mm256_insertf128_ps(_mm256_castps128_ps256(f0), f1, 0x1) };
578 }
579
580 } // namespace gmx
581
582 #endif // GMX_SIMD_IMPL_X86_AVX_256_SIMD_DOUBLE_H