Merge branch 'release-4-6' into release-5-0
[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,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.
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
58
59 /* Align a stack-based thread-local working array. */
60 static gmx_inline int *
61 prepare_table_load_buffer(const int gmx_unused *array)
62 {
63     return NULL;
64 }
65
66 #include "nbnxn_kernel_simd_utils_ref.h"
67
68 #else /* GMX_SIMD_REFERENCE */
69
70 #if defined  GMX_TARGET_X86 && !defined __MIC__
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  * 256-bit AVX use the array, but other implementations do not.
84  */
85 static gmx_inline int *
86 prepare_table_load_buffer(int gmx_unused *array)
87 {
88 #if GMX_SIMD_REAL_WIDTH >= 8 || (defined GMX_DOUBLE && GMX_SIMD_REAL_WIDTH >= 4)
89     return gmx_simd_align_i(array);
90 #else
91     return NULL;
92 #endif
93 }
94
95 #ifdef GMX_DOUBLE
96 #if GMX_SIMD_REAL_WIDTH == 2
97 #include "nbnxn_kernel_simd_utils_x86_128d.h"
98 #else
99 #include "nbnxn_kernel_simd_utils_x86_256d.h"
100 #endif
101 #else /* GMX_DOUBLE */
102 /* In single precision aligned FDV0 table loads are optimal */
103 #define TAB_FDV0
104 #if GMX_SIMD_REAL_WIDTH == 4
105 #include "nbnxn_kernel_simd_utils_x86_128s.h"
106 #else
107 #include "nbnxn_kernel_simd_utils_x86_256s.h"
108 #endif
109 #endif /* GMX_DOUBLE */
110
111 #else  /* GMX_TARGET_X86 && !__MIC__ */
112
113 #if GMX_SIMD_REAL_WIDTH > 4
114 /* For width>4 we use unaligned loads. And thus we can use the minimal stride */
115 static const int nbfp_stride = 2;
116 #else
117 static const int nbfp_stride = GMX_SIMD_REAL_WIDTH;
118 #endif
119
120 /* We use the FDV0 table layout when we can use aligned table loads */
121 #if GMX_SIMD_REAL_WIDTH == 4
122 #define TAB_FDV0
123 #endif
124
125 #ifdef GMX_SIMD_IBM_QPX
126 #include "nbnxn_kernel_simd_utils_ibm_qpx.h"
127 #endif /* GMX_SIMD_IBM_QPX */
128
129 #ifdef __MIC__
130 #include "nbnxn_kernel_simd_utils_x86_mic.h"
131 #endif
132
133 #endif /* GMX_TARGET_X86 && !__MIC__ */
134
135 #endif /* GMX_SIMD_REFERENCE */
136
137 /* If the simd width is 4, but simd4 instructions are not defined,
138  * reuse the simd real type and the four instructions we need.
139  */
140 #if GMX_SIMD_REAL_WIDTH == 4 && \
141     !((!defined GMX_DOUBLE && defined GMX_SIMD4_HAVE_FLOAT) || \
142     (defined GMX_DOUBLE && defined GMX_SIMD4_HAVE_DOUBLE))
143 #define gmx_simd4_real_t    gmx_simd_real_t
144 #define gmx_simd4_load_r    gmx_simd_load_r
145 #define gmx_simd4_store_r   gmx_simd_store_r
146 #define gmx_simd4_add_r     gmx_simd_add_r
147 #define gmx_simd4_reduce_r  gmx_simd_reduce_r
148 #endif
149
150 #ifdef UNROLLJ
151 /* Add energy register to possibly multiple terms in the energy array */
152 static gmx_inline void add_ener_grp(gmx_simd_real_t e_S, real *v, const int *offset_jj)
153 {
154     int jj;
155
156     /* We need to balance the number of store operations with
157      * the rapidly increases number of combinations of energy groups.
158      * We add to a temporary buffer for 1 i-group vs 2 j-groups.
159      */
160     for (jj = 0; jj < (UNROLLJ/2); jj++)
161     {
162         gmx_simd_real_t v_S;
163
164         v_S = gmx_simd_load_r(v+offset_jj[jj]+jj*GMX_SIMD_REAL_WIDTH);
165         gmx_simd_store_r(v+offset_jj[jj]+jj*GMX_SIMD_REAL_WIDTH, gmx_simd_add_r(v_S, e_S));
166     }
167 }
168 #endif
169
170 #if defined GMX_NBNXN_SIMD_2XNN && defined UNROLLJ
171 /* As add_ener_grp, but for two groups of UNROLLJ/2 stored in
172  * a single SIMD register.
173  */
174 static gmx_inline void
175 add_ener_grp_halves(gmx_simd_real_t e_S, real *v0, real *v1, const int *offset_jj)
176 {
177     gmx_mm_hpr e_S0, e_S1;
178     int        jj;
179
180     gmx_pr_to_2hpr(e_S, &e_S0, &e_S1);
181
182     for (jj = 0; jj < (UNROLLJ/2); jj++)
183     {
184         gmx_mm_hpr v_S;
185
186         gmx_load_hpr(&v_S, v0+offset_jj[jj]+jj*GMX_SIMD_REAL_WIDTH/2);
187         gmx_store_hpr(v0+offset_jj[jj]+jj*GMX_SIMD_REAL_WIDTH/2, gmx_add_hpr(v_S, e_S0));
188     }
189     for (jj = 0; jj < (UNROLLJ/2); jj++)
190     {
191         gmx_mm_hpr v_S;
192
193         gmx_load_hpr(&v_S, v1+offset_jj[jj]+jj*GMX_SIMD_REAL_WIDTH/2);
194         gmx_store_hpr(v1+offset_jj[jj]+jj*GMX_SIMD_REAL_WIDTH/2, gmx_add_hpr(v_S, e_S1));
195     }
196 }
197 #endif
198
199 #endif /* _nbnxn_kernel_simd_utils_h_ */