2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2013,2014, 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 _general_x86_mic_h_
37 #define _general_x86_mic_h_
39 /* This file contains the SIMD implmenetation for Intel MIC
43 #include <immintrin.h>
46 #error "Double precision isn't supported on Intel Phi yet"
49 typedef __m512 gmx_mm_ps;
50 typedef __m512 gmx_simd_real_t;
51 /* boolean SIMD register type */
52 typedef __mmask16 gmx_simd_bool_t;
53 typedef __m512i gmx_simd_int32_t;
55 #define GMX_HAVE_SIMD_MACROS
56 #define GMX_SIMD_REAL_WIDTH 16
57 #define GMX_SIMD_INT32_WIDTH 16
59 #define gmx_simd_load_r _mm512_load_ps
61 /* Set all SIMD register elements to *r */
62 static gmx_inline gmx_mm_ps
63 gmx_simd_load1_r(const real *r)
65 return _mm512_extload_ps(r, _MM_UPCONV_PS_NONE, _MM_BROADCAST_1X16, _MM_HINT_NONE);
68 #define gmx_simd_set1_r _mm512_set1_ps
69 /* Set all SIMD register elements to 0 */
70 #define gmx_simd_setzero_r _mm512_setzero_ps
71 #define gmx_simd_store_r _mm512_store_ps
73 #define gmx_simd_add_r _mm512_add_ps
74 #define gmx_simd_sub_r _mm512_sub_ps
75 #define gmx_simd_mul_r _mm512_mul_ps
77 #define GMX_SIMD_HAVE_FMA
78 #define gmx_simd_fmadd_r _mm512_fmadd_ps
79 #define gmx_simd_fnmadd_r _mm512_fnmadd_ps
81 #define gmx_simd_max_r _mm512_max_ps
83 static gmx_inline gmx_mm_ps
84 gmx_simd_blendzero_r(gmx_mm_ps a, gmx_simd_bool_t b)
86 return _mm512_mask_mov_ps(_mm512_setzero_ps(), b, a);
89 #define gmx_simd_round_r _mm512_rint_ps
91 #define GMX_SIMD_HAVE_FLOOR
92 #define gmx_simd_floor_r _mm512_floor_ps
94 /* Copy the sign of a to b, assumes b >= 0 for efficiency */
95 static gmx_inline gmx_mm_ps
96 gmx_cpsgn_nonneg_pr(gmx_mm_ps a, gmx_mm_ps b)
98 __m512 zero = _mm512_setzero_ps();
99 __m512 neg1 = _mm512_set1_ps(-1);
100 /* TODO (only bond): Bitwise operations on floating points can be done after casting to int.
101 That allows us to do it the same way as AVX which might be faster. */
102 return _mm512_mask_mul_ps(b, _mm512_cmplt_ps_mask(a, zero), b, neg1);
105 /* Very specific operation required in the non-bonded kernels */
106 static gmx_inline gmx_mm_ps
107 gmx_masknot_add_pr(gmx_simd_bool_t a, gmx_mm_ps b, gmx_mm_ps c)
109 return _mm512_mask_add_ps(b, _mm512_knot(a), b, c);
113 #define gmx_simd_cmplt_r _mm512_cmplt_ps_mask
115 /* Logical AND on SIMD booleans. */
116 #define gmx_simd_and_b _mm512_kand
118 /* Logical OR on SIMD booleans. */
119 #define gmx_simd_or_b _mm512_kor
121 /* Returns a single int (0/1) which tells if any of the booleans is True
122 It returns the full mask (not 1 for True). But given that any non-zero is True this is OK. */
123 #define gmx_simd_anytrue_b _mm512_mask2int
125 /* Conversions only used for PME table lookup */
126 static gmx_inline gmx_simd_int32_t
127 gmx_simd_cvtt_r2i(gmx_mm_ps a)
129 return _mm512_cvtfxpnt_round_adjustps_epi32(a, _MM_ROUND_MODE_DOWN, _MM_EXPADJ_NONE);
132 /* These two function only need to be approximate, Newton-Raphson iteration
133 * is used for full accuracy in gmx_simd_invsqrt_r and gmx_simd_inv_r.
135 #define gmx_simd_rsqrt_r _mm512_rsqrt23_ps
136 #define gmx_simd_rcp_r _mm512_rcp23_ps
138 #define GMX_SIMD_HAVE_EXP
139 #define gmx_simd_exp_r _mm512_exp_ps
141 #define GMX_SIMD_HAVE_ERFC
142 #define gmx_simd_erfc_r _mm512_erfc_ps
144 #define GMX_SIMD_HAVE_TRIGONOMETRIC
145 #define gmx_simd_sqrt_r _mm512_sqrt_ps
147 static gmx_inline int
148 gmx_simd_sincos_r(gmx_mm_ps a,
149 gmx_mm_ps *s, gmx_mm_ps *c)
151 /* TODO (only bond): optimize that both are calculated together.
152 Or (if if that isn't fast on MIC) don't call sincos if only one is needed. */
153 *s = _mm512_sin_ps(a);
154 *c = _mm512_cos_ps(a);
158 #define gmx_simd_acos_r _mm512_acos_ps
159 #define gmx_simd_atan2_r _mm512_atan2_ps
161 #endif /* _general_x86_mic_h_ */