5d081738f85d0359fa21ece0c35ce7eb69a8b11c
[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 #include "config.h"
41
42 /*! \brief Provides hardware-specific utility routines for the SIMD kernels.
43  *
44  * Defines all functions, typedefs, constants and macros that have
45  * explicit dependencies on the j-cluster size, precision, or SIMD
46  * width. This includes handling diagonal, Newton and topology
47  * exclusions.
48  *
49  * The functionality which depends on the j-cluster size is:
50  *   LJ-parameter lookup
51  *   force table lookup
52  *   energy group pair energy storage
53  */
54
55 #if !defined GMX_NBNXN_SIMD_2XNN && !defined GMX_NBNXN_SIMD_4XN
56 #error "Must define an NBNxN kernel flavour before including NBNxN kernel utility functions"
57 #endif
58
59 #ifdef GMX_SIMD_REFERENCE
60
61 /* Align a stack-based thread-local working array. */
62 static gmx_inline int *
63 prepare_table_load_buffer(const int gmx_unused *array)
64 {
65     return NULL;
66 }
67
68 #include "nbnxn_kernel_simd_utils_ref.h"
69
70 #else /* GMX_SIMD_REFERENCE */
71
72 #if defined  GMX_TARGET_X86 && !defined __MIC__
73 /* Include x86 SSE2 compatible SIMD functions */
74
75 /* Set the stride for the lookup of the two LJ parameters from their
76  * (padded) array. We use the minimum supported SIMD memory alignment.
77  */
78 #if defined GMX_DOUBLE
79 static const int nbfp_stride = 2;
80 #else
81 static const int nbfp_stride = 4;
82 #endif
83
84 /* Align a stack-based thread-local working array. Table loads on
85  * 256-bit AVX use the array, but other implementations do not.
86  */
87 static gmx_inline int *
88 prepare_table_load_buffer(int gmx_unused *array)
89 {
90 #if GMX_SIMD_REAL_WIDTH >= 8 || (defined GMX_DOUBLE && GMX_SIMD_REAL_WIDTH >= 4)
91     return gmx_simd_align_i(array);
92 #else
93     return NULL;
94 #endif
95 }
96
97 #ifdef GMX_DOUBLE
98 #if GMX_SIMD_REAL_WIDTH == 2
99 #include "nbnxn_kernel_simd_utils_x86_128d.h"
100 #else
101 #include "nbnxn_kernel_simd_utils_x86_256d.h"
102 #endif
103 #else /* GMX_DOUBLE */
104 /* In single precision aligned FDV0 table loads are optimal */
105 #define TAB_FDV0
106 #if GMX_SIMD_REAL_WIDTH == 4
107 #include "nbnxn_kernel_simd_utils_x86_128s.h"
108 #else
109 #include "nbnxn_kernel_simd_utils_x86_256s.h"
110 #endif
111 #endif /* GMX_DOUBLE */
112
113 #else  /* GMX_TARGET_X86 && !__MIC__ */
114
115 #if GMX_SIMD_REAL_WIDTH > 4
116 /* For width>4 we use unaligned loads. And thus we can use the minimal stride */
117 static const int nbfp_stride = 2;
118 #else
119 static const int nbfp_stride = GMX_SIMD_REAL_WIDTH;
120 #endif
121
122 /* We use the FDV0 table layout when we can use aligned table loads */
123 #if GMX_SIMD_REAL_WIDTH == 4
124 #define TAB_FDV0
125 #endif
126
127 #ifdef GMX_SIMD_IBM_QPX
128 #include "nbnxn_kernel_simd_utils_ibm_qpx.h"
129 #endif /* GMX_SIMD_IBM_QPX */
130
131 #ifdef __MIC__
132 #include "nbnxn_kernel_simd_utils_x86_mic.h"
133 #endif
134
135 #endif /* GMX_TARGET_X86 && !__MIC__ */
136
137 #endif /* GMX_SIMD_REFERENCE */
138
139 /* If the simd width is 4, but simd4 instructions are not defined,
140  * reuse the simd real type and the four instructions we need.
141  */
142 #if GMX_SIMD_REAL_WIDTH == 4 && \
143     !((!defined GMX_DOUBLE && defined GMX_SIMD4_HAVE_FLOAT) || \
144     (defined GMX_DOUBLE && defined GMX_SIMD4_HAVE_DOUBLE))
145 #define gmx_simd4_real_t    gmx_simd_real_t
146 #define gmx_simd4_load_r    gmx_simd_load_r
147 #define gmx_simd4_store_r   gmx_simd_store_r
148 #define gmx_simd4_add_r     gmx_simd_add_r
149 #define gmx_simd4_reduce_r  gmx_simd_reduce_r
150 #endif
151
152 #ifdef UNROLLJ
153 /* Add energy register to possibly multiple terms in the energy array */
154 static gmx_inline void add_ener_grp(gmx_simd_real_t e_S, real *v, const int *offset_jj)
155 {
156     int jj;
157
158     /* We need to balance the number of store operations with
159      * the rapidly increases number of combinations of energy groups.
160      * We add to a temporary buffer for 1 i-group vs 2 j-groups.
161      */
162     for (jj = 0; jj < (UNROLLJ/2); jj++)
163     {
164         gmx_simd_real_t v_S;
165
166         v_S = gmx_simd_load_r(v+offset_jj[jj]+jj*GMX_SIMD_REAL_WIDTH);
167         gmx_simd_store_r(v+offset_jj[jj]+jj*GMX_SIMD_REAL_WIDTH, gmx_simd_add_r(v_S, e_S));
168     }
169 }
170 #endif
171
172 #if defined GMX_NBNXN_SIMD_2XNN && defined UNROLLJ
173 /* As add_ener_grp, but for two groups of UNROLLJ/2 stored in
174  * a single SIMD register.
175  */
176 static gmx_inline void
177 add_ener_grp_halves(gmx_simd_real_t e_S, real *v0, real *v1, const int *offset_jj)
178 {
179     gmx_mm_hpr e_S0, e_S1;
180     int        jj;
181
182     gmx_pr_to_2hpr(e_S, &e_S0, &e_S1);
183
184     for (jj = 0; jj < (UNROLLJ/2); jj++)
185     {
186         gmx_mm_hpr v_S;
187
188         gmx_load_hpr(&v_S, v0+offset_jj[jj]+jj*GMX_SIMD_REAL_WIDTH/2);
189         gmx_store_hpr(v0+offset_jj[jj]+jj*GMX_SIMD_REAL_WIDTH/2, gmx_add_hpr(v_S, e_S0));
190     }
191     for (jj = 0; jj < (UNROLLJ/2); jj++)
192     {
193         gmx_mm_hpr v_S;
194
195         gmx_load_hpr(&v_S, v1+offset_jj[jj]+jj*GMX_SIMD_REAL_WIDTH/2);
196         gmx_store_hpr(v1+offset_jj[jj]+jj*GMX_SIMD_REAL_WIDTH/2, gmx_add_hpr(v_S, e_S1));
197     }
198 }
199 #endif
200
201 #endif /* _nbnxn_kernel_simd_utils_h_ */