9fde8f06e7c5caa2f9c53dbf7e4dd8db5c18d40d
[alexxy/gromacs.git] / src / mdlib / nbnxn_kernels / nbnxn_kernel_simd_utils.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2012, The GROMACS Development Team
6  * Copyright (c) 2012,2013, by the GROMACS development team, led by
7  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
8  * others, as listed in the AUTHORS file in the top-level source
9  * directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 #ifndef _nbnxn_kernel_simd_utils_h_
38 #define _nbnxn_kernel_simd_utils_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 /* Set the stride for the lookup of the two LJ parameters from their
60  * (padded) array.
61  * Note that currently only arrays with stride 2 and 4 are available.
62  * Since the reference code does not require alignment, we can always use 2.
63  */
64 static const int nbfp_stride = 2;
65
66 /* Align a stack-based thread-local working array. */
67 static gmx_inline int *
68 prepare_table_load_buffer(const int *array)
69 {
70     return NULL;
71 }
72
73 #include "nbnxn_kernel_simd_utils_ref.h"
74
75 #else /* GMX_SIMD_REFERENCE_PLAIN_C */
76
77 #ifdef GMX_X86_SSE2
78 /* Include x86 SSE2 compatible SIMD functions */
79
80 /* Set the stride for the lookup of the two LJ parameters from their
81  * (padded) array. We use the minimum supported SIMD memory alignment.
82  */
83 #if defined GMX_DOUBLE
84 static const int nbfp_stride = 2;
85 #else
86 static const int nbfp_stride = 4;
87 #endif
88
89 /* Align a stack-based thread-local working array. Table loads on
90  * full-width AVX_256 use the array, but other implementations do
91  * not. */
92 static gmx_inline int *
93 prepare_table_load_buffer(const int *array)
94 {
95 #if defined GMX_X86_AVX_256 && !defined GMX_USE_HALF_WIDTH_SIMD_HERE
96     return gmx_simd_align_int(array);
97 #else
98     return NULL;
99 #endif
100 }
101
102 #if defined GMX_X86_AVX_256 && !defined GMX_USE_HALF_WIDTH_SIMD_HERE
103
104 /* With full AVX-256 SIMD, half SIMD-width table loads are optimal */
105 #if GMX_SIMD_WIDTH_HERE == 8
106 #define TAB_FDV0
107 #endif
108
109 /*
110 Berk, 2xnn.c had the following code, but I think it is safe to remove now, given the code immediately above.
111
112 #if defined GMX_X86_AVX_256 && !defined GMX_DOUBLE
113 / * AVX-256 single precision 2x(4+4) kernel,
114  * we can do half SIMD-width aligned FDV0 table loads.
115  * /
116 #define TAB_FDV0
117 #endif
118 */
119
120 #ifdef GMX_DOUBLE
121 #include "nbnxn_kernel_simd_utils_x86_256d.h"
122 #else  /* GMX_DOUBLE */
123 #include "nbnxn_kernel_simd_utils_x86_256s.h"
124 #endif /* GMX_DOUBLE */
125
126 #else  /* defined GMX_X86_AVX_256 && !defined GMX_USE_HALF_WIDTH_SIMD_HERE */
127
128 /* We use the FDV0 table layout when we can use aligned table loads */
129 #if GMX_SIMD_WIDTH_HERE == 4
130 #define TAB_FDV0
131 #endif
132
133 #ifdef GMX_DOUBLE
134 #include "nbnxn_kernel_simd_utils_x86_128d.h"
135 #else  /* GMX_DOUBLE */
136 #include "nbnxn_kernel_simd_utils_x86_128s.h"
137 #endif /* GMX_DOUBLE */
138
139 #endif /* defined GMX_X86_AVX_256 && !defined GMX_USE_HALF_WIDTH_SIMD_HERE */
140
141 #else  /* GMX_X86_SSE2 */
142
143 #if GMX_SIMD_WIDTH_HERE > 4
144 static const int nbfp_stride = 4;
145 #else
146 static const int nbfp_stride = GMX_SIMD_WIDTH_HERE;
147 #endif
148
149 /* We use the FDV0 table layout when we can use aligned table loads */
150 #if GMX_SIMD_WIDTH_HERE == 4
151 #define TAB_FDV0
152 #endif
153
154 #ifdef GMX_CPU_ACCELERATION_IBM_QPX
155 #include "nbnxn_kernel_simd_utils_ibm_qpx.h"
156 #endif /* GMX_CPU_ACCELERATION_IBM_QPX */
157
158 #endif /* GMX_X86_SSE2 */
159 #endif /* GMX_SIMD_REFERENCE_PLAIN_C */
160
161
162 #ifdef UNROLLJ
163 /* Add energy register to possibly multiple terms in the energy array */
164 static inline void add_ener_grp(gmx_mm_pr e_S, real *v, const int *offset_jj)
165 {
166     int jj;
167
168     /* We need to balance the number of store operations with
169      * the rapidly increases number of combinations of energy groups.
170      * We add to a temporary buffer for 1 i-group vs 2 j-groups.
171      */
172     for (jj = 0; jj < (UNROLLJ/2); jj++)
173     {
174         gmx_mm_pr v_S;
175
176         v_S = gmx_load_pr(v+offset_jj[jj]+jj*GMX_SIMD_WIDTH_HERE);
177         gmx_store_pr(v+offset_jj[jj]+jj*GMX_SIMD_WIDTH_HERE, gmx_add_pr(v_S, e_S));
178     }
179 }
180 #endif
181
182 #if defined GMX_NBNXN_SIMD_2XNN && defined UNROLLJ
183 /* As add_ener_grp, but for two groups of UNROLLJ/2 stored in
184  * a single SIMD register.
185  */
186 static inline void
187 add_ener_grp_halves(gmx_mm_pr e_S, real *v0, real *v1, const int *offset_jj)
188 {
189     gmx_mm_hpr e_S0, e_S1;
190     int        jj;
191
192     gmx_pr_to_2hpr(e_S, &e_S0, &e_S1);
193
194     for (jj = 0; jj < (UNROLLJ/2); jj++)
195     {
196         gmx_mm_hpr v_S;
197
198         gmx_load_hpr(&v_S, v0+offset_jj[jj]+jj*GMX_SIMD_WIDTH_HERE/2);
199         gmx_store_hpr(v0+offset_jj[jj]+jj*GMX_SIMD_WIDTH_HERE/2, gmx_add_hpr(v_S, e_S0));
200     }
201     for (jj = 0; jj < (UNROLLJ/2); jj++)
202     {
203         gmx_mm_hpr v_S;
204
205         gmx_load_hpr(&v_S, v1+offset_jj[jj]+jj*GMX_SIMD_WIDTH_HERE/2);
206         gmx_store_hpr(v1+offset_jj[jj]+jj*GMX_SIMD_WIDTH_HERE/2, gmx_add_hpr(v_S, e_S1));
207     }
208 }
209 #endif
210
211 #endif /* _nbnxn_kernel_simd_utils_h_ */