Apply clang-format to source tree
[alexxy/gromacs.git] / src / gromacs / simd / impl_x86_avx_128_fma / impl_x86_avx_128_fma_util_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
36 #ifndef GMX_SIMD_IMPL_X86_AVX_128_FMA_UTIL_FLOAT_H
37 #define GMX_SIMD_IMPL_X86_AVX_128_FMA_UTIL_FLOAT_H
38
39 #include "config.h"
40
41 #include <cassert>
42 #include <cstddef>
43 #include <cstdint>
44
45 #include <immintrin.h>
46 #include <x86intrin.h>
47
48 #include "gromacs/simd/impl_x86_sse4_1/impl_x86_sse4_1_util_float.h"
49
50 namespace gmx
51 {
52
53 /* In the old group kernels we found it more efficient to transpose the data to store rather
54  * than using maskload and maskstore. It might be worth to test again, but for now we assume
55  * this is still the case, and rely on those version inherited from the SSE2 header.
56  *
57  * It is also worth testing if changing _mm_shuffle_ps() to _mm_permute_ps() could improve
58  * throughput just-so-slightly.
59  */
60
61 static inline void gmx_simdcall expandScalarsToTriplets(SimdFloat  scalar,
62                                                         SimdFloat* triplets0,
63                                                         SimdFloat* triplets1,
64                                                         SimdFloat* triplets2)
65 {
66     triplets0->simdInternal_ = _mm_permute_ps(scalar.simdInternal_, _MM_SHUFFLE(1, 0, 0, 0));
67     triplets1->simdInternal_ = _mm_permute_ps(scalar.simdInternal_, _MM_SHUFFLE(2, 2, 1, 1));
68     triplets2->simdInternal_ = _mm_permute_ps(scalar.simdInternal_, _MM_SHUFFLE(3, 3, 3, 2));
69 }
70
71 static inline float gmx_simdcall reduceIncr4ReturnSum(float* m, SimdFloat v0, SimdFloat v1, SimdFloat v2, SimdFloat v3)
72 {
73     _MM_TRANSPOSE4_PS(v0.simdInternal_, v1.simdInternal_, v2.simdInternal_, v3.simdInternal_);
74     v0.simdInternal_ = _mm_add_ps(v0.simdInternal_, v1.simdInternal_);
75     v2.simdInternal_ = _mm_add_ps(v2.simdInternal_, v3.simdInternal_);
76     v0.simdInternal_ = _mm_add_ps(v0.simdInternal_, v2.simdInternal_);
77
78     assert(std::size_t(m) % 16 == 0);
79
80     v2.simdInternal_ = _mm_add_ps(v0.simdInternal_, _mm_load_ps(m));
81     _mm_store_ps(m, v2.simdInternal_);
82
83     __m128 b = _mm_add_ps(v0.simdInternal_, _mm_permute_ps(v0.simdInternal_, _MM_SHUFFLE(1, 0, 3, 2)));
84     b        = _mm_add_ss(b, _mm_permute_ps(b, _MM_SHUFFLE(0, 3, 2, 1)));
85     return *reinterpret_cast<float*>(&b);
86 }
87
88 } // namespace gmx
89
90 #endif // GMX_SIMD_IMPL_X86_AVX_128_FMA_UTIL_FLOAT_H