2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,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.
35 #ifndef _nbnxn_kernel_simd_utils_h_
36 #define _nbnxn_kernel_simd_utils_h_
38 #include "gromacs/legacyheaders/types/simple.h"
40 /*! \brief Provides hardware-specific utility routines for the SIMD kernels.
42 * Defines all functions, typedefs, constants and macros that have
43 * explicit dependencies on the j-cluster size, precision, or SIMD
44 * width. This includes handling diagonal, Newton and topology
47 * The functionality which depends on the j-cluster size is:
50 * energy group pair energy storage
53 #if !defined GMX_NBNXN_SIMD_2XNN && !defined GMX_NBNXN_SIMD_4XN
54 #error "Must define an NBNxN kernel flavour before including NBNxN kernel utility functions"
57 #ifdef GMX_SIMD_REFERENCE
59 /* Align a stack-based thread-local working array. */
60 static gmx_inline int *
61 prepare_table_load_buffer(const int gmx_unused *array)
66 #include "nbnxn_kernel_simd_utils_ref.h"
68 #else /* GMX_SIMD_REFERENCE */
70 #if defined GMX_SIMD_X86_SSE2_OR_HIGHER && !defined __MIC__
71 /* Include x86 SSE2 compatible SIMD functions */
73 /* Set the stride for the lookup of the two LJ parameters from their
74 * (padded) array. We use the minimum supported SIMD memory alignment.
76 #if defined GMX_DOUBLE
77 static const int nbfp_stride = 2;
79 static const int nbfp_stride = 4;
82 /* Align a stack-based thread-local working array. Table loads on
83 * full-width AVX_256 use the array, but other implementations do
85 static gmx_inline int *
86 prepare_table_load_buffer(const int gmx_unused *array)
88 #if defined GMX_SIMD_X86_AVX_256_OR_HIGHER && !defined GMX_USE_HALF_WIDTH_SIMD_HERE
89 return gmx_simd_align_i(array);
95 #if defined GMX_SIMD_X86_AVX_256_OR_HIGHER && !defined GMX_USE_HALF_WIDTH_SIMD_HERE
97 /* With full AVX-256 SIMD, half SIMD-width table loads are optimal */
98 #if GMX_SIMD_REAL_WIDTH == 8
102 #include "nbnxn_kernel_simd_utils_x86_256d.h"
103 #else /* GMX_DOUBLE */
104 #include "nbnxn_kernel_simd_utils_x86_256s.h"
105 #endif /* GMX_DOUBLE */
107 #else /* defined GMX_SIMD_X86_AVX_256_OR_HIGHER && !defined GMX_USE_HALF_WIDTH_SIMD_HERE */
109 /* We use the FDV0 table layout when we can use aligned table loads */
110 #if GMX_SIMD_REAL_WIDTH == 4
115 #include "nbnxn_kernel_simd_utils_x86_128d.h"
116 #else /* GMX_DOUBLE */
117 #include "nbnxn_kernel_simd_utils_x86_128s.h"
118 #endif /* GMX_DOUBLE */
120 #endif /* defined GMX_SIMD_X86_AVX_256_OR_HIGHER && !defined GMX_USE_HALF_WIDTH_SIMD_HERE */
122 #else /* GMX_SIMD_X86_SSE2_OR_HIGHER */
124 #if GMX_SIMD_REAL_WIDTH > 4
125 /* For width>4 we use unaligned loads. And thus we can use the minimal stride */
126 static const int nbfp_stride = 2;
128 static const int nbfp_stride = GMX_SIMD_REAL_WIDTH;
131 /* We use the FDV0 table layout when we can use aligned table loads */
132 #if GMX_SIMD_REAL_WIDTH == 4
136 #ifdef GMX_SIMD_IBM_QPX
137 #include "nbnxn_kernel_simd_utils_ibm_qpx.h"
138 #endif /* GMX_SIMD_IBM_QPX */
141 #include "nbnxn_kernel_simd_utils_x86_mic.h"
144 #endif /* GMX_SIMD_X86_SSE2_OR_HIGHER */
145 #endif /* GMX_SIMD_REFERENCE */
147 #if GMX_SIMD_REAL_WIDTH == 4
148 #define gmx_mm_pr4 gmx_simd_real_t
149 #define gmx_load_pr4 gmx_simd_load_r
150 #define gmx_store_pr4 gmx_simd_store_r
151 #define gmx_add_pr4 gmx_simd_add_r
154 #ifndef HAVE_GMX_SUM_SIMD /* should be defined for arch with hardware reduce */
155 static gmx_inline real
156 gmx_sum_simd2(gmx_simd_real_t x, real* b)
158 gmx_simd_store_r(b, x);
162 #if GMX_SIMD_REAL_WIDTH >= 4
163 static gmx_inline real
164 gmx_sum_simd4(gmx_mm_pr4 x, real* b)
167 return b[0]+b[1]+b[2]+b[3];
171 #if GMX_SIMD_REAL_WIDTH == 2
172 static gmx_inline real gmx_sum_simd(gmx_simd_real_t x, real* b)
174 gmx_simd_store_r(b, x);
177 #elif GMX_SIMD_REAL_WIDTH == 4
178 static gmx_inline real gmx_sum_simd(gmx_simd_real_t x, real* b)
180 gmx_simd_store_r(b, x);
181 return b[0]+b[1]+b[2]+b[3];
183 #elif GMX_SIMD_REAL_WIDTH == 8
184 static gmx_inline real gmx_sum_simd(gmx_simd_real_t x, real* b)
186 gmx_simd_store_r(b, x);
187 return b[0]+b[1]+b[2]+b[3]+b[4]+b[5]+b[6]+b[7];
189 #elif GMX_SIMD_REAL_WIDTH == 16
190 /* This is getting ridiculous, SIMD horizontal adds would help,
191 * but this is not performance critical (only used to reduce energies)
193 static gmx_inline real gmx_sum_simd(gmx_simd_real_t x, real* b)
195 gmx_simd_store_r(b, x);
196 return b[0]+b[1]+b[2]+b[3]+b[4]+b[5]+b[6]+b[7]+b[8]+b[9]+b[10]+b[11]+b[12]+b[13]+b[14]+b[15];
199 #error "unsupported kernel configuration"
201 #endif //HAVE_GMX_SUM_SIMD
204 /* Add energy register to possibly multiple terms in the energy array */
205 static inline void add_ener_grp(gmx_simd_real_t e_S, real *v, const int *offset_jj)
209 /* We need to balance the number of store operations with
210 * the rapidly increases number of combinations of energy groups.
211 * We add to a temporary buffer for 1 i-group vs 2 j-groups.
213 for (jj = 0; jj < (UNROLLJ/2); jj++)
217 v_S = gmx_simd_load_r(v+offset_jj[jj]+jj*GMX_SIMD_REAL_WIDTH);
218 gmx_simd_store_r(v+offset_jj[jj]+jj*GMX_SIMD_REAL_WIDTH, gmx_simd_add_r(v_S, e_S));
223 #if defined GMX_NBNXN_SIMD_2XNN && defined UNROLLJ
224 /* As add_ener_grp, but for two groups of UNROLLJ/2 stored in
225 * a single SIMD register.
228 add_ener_grp_halves(gmx_simd_real_t e_S, real *v0, real *v1, const int *offset_jj)
230 gmx_mm_hpr e_S0, e_S1;
233 gmx_pr_to_2hpr(e_S, &e_S0, &e_S1);
235 for (jj = 0; jj < (UNROLLJ/2); jj++)
239 gmx_load_hpr(&v_S, v0+offset_jj[jj]+jj*GMX_SIMD_REAL_WIDTH/2);
240 gmx_store_hpr(v0+offset_jj[jj]+jj*GMX_SIMD_REAL_WIDTH/2, gmx_add_hpr(v_S, e_S0));
242 for (jj = 0; jj < (UNROLLJ/2); jj++)
246 gmx_load_hpr(&v_S, v1+offset_jj[jj]+jj*GMX_SIMD_REAL_WIDTH/2);
247 gmx_store_hpr(v1+offset_jj[jj]+jj*GMX_SIMD_REAL_WIDTH/2, gmx_add_hpr(v_S, e_S1));
252 #endif /* _nbnxn_kernel_simd_utils_h_ */