Replace SUM_SIMD with gmx_sum_simd
[alexxy/gromacs.git] / src / gromacs / mdlib / nbnxn_kernels / nbnxn_kernel_simd_utils.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013, 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 #ifndef _nbnxn_kernel_simd_utils_h_
36 #define _nbnxn_kernel_simd_utils_h_
37
38 #include "gromacs/legacyheaders/types/simple.h"
39
40 /*! \brief Provides hardware-specific utility routines for the SIMD kernels.
41  *
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
45  * exclusions.
46  *
47  * The functionality which depends on the j-cluster size is:
48  *   LJ-parameter lookup
49  *   force table lookup
50  *   energy group pair energy storage
51  */
52
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"
55 #endif
56
57 #ifdef GMX_SIMD_REFERENCE_PLAIN_C
58
59 /* Align a stack-based thread-local working array. */
60 static gmx_inline int *
61 prepare_table_load_buffer(const int *array)
62 {
63     return NULL;
64 }
65
66 #include "nbnxn_kernel_simd_utils_ref.h"
67
68 #else /* GMX_SIMD_REFERENCE_PLAIN_C */
69
70 #ifdef GMX_X86_SSE2
71 /* Include x86 SSE2 compatible SIMD functions */
72
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.
75  */
76 #if defined GMX_DOUBLE
77 static const int nbfp_stride = 2;
78 #else
79 static const int nbfp_stride = 4;
80 #endif
81
82 /* Align a stack-based thread-local working array. Table loads on
83  * full-width AVX_256 use the array, but other implementations do
84  * not. */
85 static gmx_inline int *
86 prepare_table_load_buffer(const int gmx_unused *array)
87 {
88 #if defined GMX_X86_AVX_256 && !defined GMX_USE_HALF_WIDTH_SIMD_HERE
89     return gmx_simd_align_int(array);
90 #else
91     return NULL;
92 #endif
93 }
94
95 #if defined GMX_X86_AVX_256 && !defined GMX_USE_HALF_WIDTH_SIMD_HERE
96
97 /* With full AVX-256 SIMD, half SIMD-width table loads are optimal */
98 #if GMX_SIMD_WIDTH_HERE == 8
99 #define TAB_FDV0
100 #endif
101
102 #ifdef GMX_DOUBLE
103 #include "nbnxn_kernel_simd_utils_x86_256d.h"
104 #else  /* GMX_DOUBLE */
105 #include "nbnxn_kernel_simd_utils_x86_256s.h"
106 #endif /* GMX_DOUBLE */
107
108 #else  /* defined GMX_X86_AVX_256 && !defined GMX_USE_HALF_WIDTH_SIMD_HERE */
109
110 /* We use the FDV0 table layout when we can use aligned table loads */
111 #if GMX_SIMD_WIDTH_HERE == 4
112 #define TAB_FDV0
113 #endif
114
115 #ifdef GMX_DOUBLE
116 #include "nbnxn_kernel_simd_utils_x86_128d.h"
117 #else  /* GMX_DOUBLE */
118 #include "nbnxn_kernel_simd_utils_x86_128s.h"
119 #endif /* GMX_DOUBLE */
120
121 #endif /* defined GMX_X86_AVX_256 && !defined GMX_USE_HALF_WIDTH_SIMD_HERE */
122
123 #else  /* GMX_X86_SSE2 */
124
125 #if GMX_SIMD_WIDTH_HERE > 4
126 static const int nbfp_stride = 4;
127 #else
128 static const int nbfp_stride = GMX_SIMD_WIDTH_HERE;
129 #endif
130
131 /* We use the FDV0 table layout when we can use aligned table loads */
132 #if GMX_SIMD_WIDTH_HERE == 4
133 #define TAB_FDV0
134 #endif
135
136 #ifdef GMX_CPU_ACCELERATION_IBM_QPX
137 #include "nbnxn_kernel_simd_utils_ibm_qpx.h"
138 #endif /* GMX_CPU_ACCELERATION_IBM_QPX */
139
140 #endif /* GMX_X86_SSE2 */
141 #endif /* GMX_SIMD_REFERENCE_PLAIN_C */
142
143 #if GMX_SIMD_WIDTH_HERE == 4 && !defined GMX_SIMD_REFERENCE_PLAIN_C
144 #define gmx_mm_pr4    gmx_mm_pr
145 #define gmx_load_pr4  gmx_load_pr
146 #define gmx_store_pr4 gmx_store_pr
147 #define gmx_add_pr4   gmx_add_pr
148 #endif
149
150 #ifndef HAVE_GMX_SUM_SIMD /* should be defined for arch with hardware reduce */
151 static gmx_inline real
152 gmx_sum_simd2(gmx_mm_pr x, real* b)
153 {
154     gmx_store_pr(b, x);
155     return b[0]+b[1];
156 }
157
158 #if GMX_SIMD_WIDTH_HERE>=4
159 static gmx_inline real
160 gmx_sum_simd4(gmx_mm_pr4 x, real* b)
161 {
162     gmx_store_pr4(b, x);
163     return b[0]+b[1]+b[2]+b[3];
164 }
165 #endif
166
167 #if GMX_SIMD_WIDTH_HERE == 2
168 static gmx_inline real gmx_sum_simd(gmx_mm_pr x, real* b)
169 {
170     gmx_store_pr(b, x);
171     return b[0]+b[1];
172 }
173 #elif GMX_SIMD_WIDTH_HERE == 4
174 static gmx_inline real gmx_sum_simd(gmx_mm_pr x, real* b)
175 {
176     gmx_store_pr(b, x);
177     return b[0]+b[1]+b[2]+b[3];
178 }
179 #elif GMX_SIMD_WIDTH_HERE == 8
180 static gmx_inline real gmx_sum_simd(gmx_mm_pr x, real* b)
181 {
182     gmx_store_pr(b, x);
183     return b[0]+b[1]+b[2]+b[3]+b[4]+b[5]+b[6]+b[7];
184 }
185 #elif GMX_SIMD_WIDTH_HERE == 16
186 /* This is getting ridiculous, SIMD horizontal adds would help,
187  * but this is not performance critical (only used to reduce energies)
188  */
189 static gmx_inline real gmx_sum_simd(gmx_mm_pr x, real* b)
190 {
191     gmx_store_pr(b, x);
192     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];
193 }
194 #else
195 #error "unsupported kernel configuration"
196 #endif
197 #endif //HAVE_GMX_SUM_SIMD
198
199 #ifdef UNROLLJ
200 /* Add energy register to possibly multiple terms in the energy array */
201 static inline void add_ener_grp(gmx_mm_pr e_S, real *v, const int *offset_jj)
202 {
203     int jj;
204
205     /* We need to balance the number of store operations with
206      * the rapidly increases number of combinations of energy groups.
207      * We add to a temporary buffer for 1 i-group vs 2 j-groups.
208      */
209     for (jj = 0; jj < (UNROLLJ/2); jj++)
210     {
211         gmx_mm_pr v_S;
212
213         v_S = gmx_load_pr(v+offset_jj[jj]+jj*GMX_SIMD_WIDTH_HERE);
214         gmx_store_pr(v+offset_jj[jj]+jj*GMX_SIMD_WIDTH_HERE, gmx_add_pr(v_S, e_S));
215     }
216 }
217 #endif
218
219 #if defined GMX_NBNXN_SIMD_2XNN && defined UNROLLJ
220 /* As add_ener_grp, but for two groups of UNROLLJ/2 stored in
221  * a single SIMD register.
222  */
223 static inline void
224 add_ener_grp_halves(gmx_mm_pr e_S, real *v0, real *v1, const int *offset_jj)
225 {
226     gmx_mm_hpr e_S0, e_S1;
227     int        jj;
228
229     gmx_pr_to_2hpr(e_S, &e_S0, &e_S1);
230
231     for (jj = 0; jj < (UNROLLJ/2); jj++)
232     {
233         gmx_mm_hpr v_S;
234
235         gmx_load_hpr(&v_S, v0+offset_jj[jj]+jj*GMX_SIMD_WIDTH_HERE/2);
236         gmx_store_hpr(v0+offset_jj[jj]+jj*GMX_SIMD_WIDTH_HERE/2, gmx_add_hpr(v_S, e_S0));
237     }
238     for (jj = 0; jj < (UNROLLJ/2); jj++)
239     {
240         gmx_mm_hpr v_S;
241
242         gmx_load_hpr(&v_S, v1+offset_jj[jj]+jj*GMX_SIMD_WIDTH_HERE/2);
243         gmx_store_hpr(v1+offset_jj[jj]+jj*GMX_SIMD_WIDTH_HERE/2, gmx_add_hpr(v_S, e_S1));
244     }
245 }
246 #endif
247
248 #endif /* _nbnxn_kernel_simd_utils_h_ */