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 #include "nbnxn_kernel_simd_utils_ref.h"
60
61 #else /* GMX_SIMD_REFERENCE_PLAIN_C */
62
63 #ifdef GMX_X86_SSE2
64 /* Include x86 SSE2 compatible SIMD functions */
65
66 /* Set the stride for the lookup of the two LJ parameters from their
67    (padded) array. Only strides of 2 and 4 are currently supported. */
68 #if defined GMX_NBNXN_SIMD_2XNN
69 static const int nbfp_stride = 4;
70 #elif defined GMX_DOUBLE
71 static const int nbfp_stride = 2;
72 #else
73 static const int nbfp_stride = 4;
74 #endif
75
76 /* Align a stack-based thread-local working array. Table loads on
77  * full-width AVX_256 use the array, but other implementations do
78  * not. */
79 static gmx_inline int *
80 prepare_table_load_buffer(const int *array)
81 {
82 #if defined GMX_X86_AVX_256 && !defined GMX_USE_HALF_WIDTH_SIMD_HERE
83     return gmx_simd_align_int(array);
84 #else
85     return NULL;
86 #endif
87 }
88
89 #if defined GMX_X86_AVX_256 && !defined GMX_USE_HALF_WIDTH_SIMD_HERE
90
91 /* With full AVX-256 SIMD, half SIMD-width table loads are optimal */
92 #if GMX_SIMD_WIDTH_HERE == 8
93 #define TAB_FDV0
94 #endif
95
96 /*
97 Berk, 2xnn.c had the following code, but I think it is safe to remove now, given the code immediately above.
98
99 #if defined GMX_X86_AVX_256 && !defined GMX_DOUBLE
100 / * AVX-256 single precision 2x(4+4) kernel,
101  * we can do half SIMD-width aligned FDV0 table loads.
102  * /
103 #define TAB_FDV0
104 #endif
105 */
106
107 #ifdef GMX_DOUBLE
108 #include "nbnxn_kernel_simd_utils_x86_256d.h"
109 #else  /* GMX_DOUBLE */
110 #include "nbnxn_kernel_simd_utils_x86_256s.h"
111 #endif /* GMX_DOUBLE */
112
113 #else  /* defined GMX_X86_AVX_256 && !defined GMX_USE_HALF_WIDTH_SIMD_HERE */
114
115 /* We use the FDV0 table layout when we can use aligned table loads */
116 #if GMX_SIMD_WIDTH_HERE == 4
117 #define TAB_FDV0
118 #endif
119
120 #ifdef GMX_DOUBLE
121 #include "nbnxn_kernel_simd_utils_x86_128d.h"
122 #else  /* GMX_DOUBLE */
123 #include "nbnxn_kernel_simd_utils_x86_128s.h"
124 #endif /* GMX_DOUBLE */
125
126 #endif /* defined GMX_X86_AVX_256 && !defined GMX_USE_HALF_WIDTH_SIMD_HERE */
127
128 #else  /* GMX_X86_SSE2 */
129
130 #if GMX_SIMD_WIDTH_HERE > 4
131 static const int nbfp_stride = 4;
132 #else
133 static const int nbfp_stride = GMX_SIMD_WIDTH_HERE;
134 #endif
135
136 #endif /* GMX_X86_SSE2 */
137 #endif /* GMX_SIMD_REFERENCE_PLAIN_C */
138
139
140 #ifdef UNROLLJ
141 /* Add energy register to possibly multiple terms in the energy array */
142 static inline void add_ener_grp(gmx_mm_pr e_S, real *v, const int *offset_jj)
143 {
144     int jj;
145
146     /* We need to balance the number of store operations with
147      * the rapidly increases number of combinations of energy groups.
148      * We add to a temporary buffer for 1 i-group vs 2 j-groups.
149      */
150     for (jj = 0; jj < (UNROLLJ/2); jj++)
151     {
152         gmx_mm_pr v_S;
153
154         v_S = gmx_load_pr(v+offset_jj[jj]+jj*GMX_SIMD_WIDTH_HERE);
155         gmx_store_pr(v+offset_jj[jj]+jj*GMX_SIMD_WIDTH_HERE, gmx_add_pr(v_S, e_S));
156     }
157 }
158 #endif
159
160 #if defined GMX_NBNXN_SIMD_2XNN && defined UNROLLJ
161 /* As add_ener_grp, but for two groups of UNROLLJ/2 stored in
162  * a single SIMD register.
163  */
164 static inline void
165 add_ener_grp_halves(gmx_mm_pr e_S, real *v0, real *v1, const int *offset_jj)
166 {
167     gmx_mm_hpr e_S0, e_S1;
168     int        jj;
169
170     gmx_pr_to_2hpr(e_S, &e_S0, &e_S1);
171
172     for (jj = 0; jj < (UNROLLJ/2); jj++)
173     {
174         gmx_mm_hpr v_S;
175
176         gmx_load_hpr(&v_S, v0+offset_jj[jj]+jj*GMX_SIMD_WIDTH_HERE/2);
177         gmx_store_hpr(v0+offset_jj[jj]+jj*GMX_SIMD_WIDTH_HERE/2, gmx_add_hpr(v_S, e_S0));
178     }
179     for (jj = 0; jj < (UNROLLJ/2); jj++)
180     {
181         gmx_mm_hpr v_S;
182
183         gmx_load_hpr(&v_S, v1+offset_jj[jj]+jj*GMX_SIMD_WIDTH_HERE/2);
184         gmx_store_hpr(v1+offset_jj[jj]+jj*GMX_SIMD_WIDTH_HERE/2, gmx_add_hpr(v_S, e_S1));
185     }
186 }
187 #endif
188
189 #endif /* _nbnxn_kernel_simd_utils_h_ */