Merge branch release-4-6
[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) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2012, The GROMACS Development Team
6  * Copyright (c) 2012, 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 /* 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 *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 /*
103 Berk, 2xnn.c had the following code, but I think it is safe to remove now, given the code immediately above.
104
105 #if defined GMX_X86_AVX_256 && !defined GMX_DOUBLE
106 / * AVX-256 single precision 2x(4+4) kernel,
107  * we can do half SIMD-width aligned FDV0 table loads.
108  * /
109 #define TAB_FDV0
110 #endif
111 */
112
113 #ifdef GMX_DOUBLE
114 #include "nbnxn_kernel_simd_utils_x86_256d.h"
115 #else  /* GMX_DOUBLE */
116 #include "nbnxn_kernel_simd_utils_x86_256s.h"
117 #endif /* GMX_DOUBLE */
118
119 #else  /* defined GMX_X86_AVX_256 && !defined GMX_USE_HALF_WIDTH_SIMD_HERE */
120
121 /* We use the FDV0 table layout when we can use aligned table loads */
122 #if GMX_SIMD_WIDTH_HERE == 4
123 #define TAB_FDV0
124 #endif
125
126 #ifdef GMX_DOUBLE
127 #include "nbnxn_kernel_simd_utils_x86_128d.h"
128 #else  /* GMX_DOUBLE */
129 #include "nbnxn_kernel_simd_utils_x86_128s.h"
130 #endif /* GMX_DOUBLE */
131
132 #endif /* defined GMX_X86_AVX_256 && !defined GMX_USE_HALF_WIDTH_SIMD_HERE */
133
134 #else  /* GMX_X86_SSE2 */
135
136 #if GMX_SIMD_WIDTH_HERE > 4
137 static const int nbfp_stride = 4;
138 #else
139 static const int nbfp_stride = GMX_SIMD_WIDTH_HERE;
140 #endif
141
142 /* We use the FDV0 table layout when we can use aligned table loads */
143 #if GMX_SIMD_WIDTH_HERE == 4
144 #define TAB_FDV0
145 #endif
146
147 #ifdef GMX_CPU_ACCELERATION_IBM_QPX
148 #include "nbnxn_kernel_simd_utils_ibm_qpx.h"
149 #endif /* GMX_CPU_ACCELERATION_IBM_QPX */
150
151 #endif /* GMX_X86_SSE2 */
152 #endif /* GMX_SIMD_REFERENCE_PLAIN_C */
153
154
155 #ifdef UNROLLJ
156 /* Add energy register to possibly multiple terms in the energy array */
157 static inline void add_ener_grp(gmx_mm_pr e_S, real *v, const int *offset_jj)
158 {
159     int jj;
160
161     /* We need to balance the number of store operations with
162      * the rapidly increases number of combinations of energy groups.
163      * We add to a temporary buffer for 1 i-group vs 2 j-groups.
164      */
165     for (jj = 0; jj < (UNROLLJ/2); jj++)
166     {
167         gmx_mm_pr v_S;
168
169         v_S = gmx_load_pr(v+offset_jj[jj]+jj*GMX_SIMD_WIDTH_HERE);
170         gmx_store_pr(v+offset_jj[jj]+jj*GMX_SIMD_WIDTH_HERE, gmx_add_pr(v_S, e_S));
171     }
172 }
173 #endif
174
175 #if defined GMX_NBNXN_SIMD_2XNN && defined UNROLLJ
176 /* As add_ener_grp, but for two groups of UNROLLJ/2 stored in
177  * a single SIMD register.
178  */
179 static inline void
180 add_ener_grp_halves(gmx_mm_pr e_S, real *v0, real *v1, const int *offset_jj)
181 {
182     gmx_mm_hpr e_S0, e_S1;
183     int        jj;
184
185     gmx_pr_to_2hpr(e_S, &e_S0, &e_S1);
186
187     for (jj = 0; jj < (UNROLLJ/2); jj++)
188     {
189         gmx_mm_hpr v_S;
190
191         gmx_load_hpr(&v_S, v0+offset_jj[jj]+jj*GMX_SIMD_WIDTH_HERE/2);
192         gmx_store_hpr(v0+offset_jj[jj]+jj*GMX_SIMD_WIDTH_HERE/2, gmx_add_hpr(v_S, e_S0));
193     }
194     for (jj = 0; jj < (UNROLLJ/2); jj++)
195     {
196         gmx_mm_hpr v_S;
197
198         gmx_load_hpr(&v_S, v1+offset_jj[jj]+jj*GMX_SIMD_WIDTH_HERE/2);
199         gmx_store_hpr(v1+offset_jj[jj]+jj*GMX_SIMD_WIDTH_HERE/2, gmx_add_hpr(v_S, e_S1));
200     }
201 }
202 #endif
203
204 #endif /* _nbnxn_kernel_simd_utils_h_ */