Updated dlist.c to recognize more atom names.
[alexxy/gromacs.git] / src / mdlib / nbnxn_kernels / nbnxn_kernel_simd_utils_x86_128d.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_x86_128d_h_
38 #define _nbnxn_kernel_simd_utils_x86_128d_h_
39
40 /* This files contains all functions/macros for the SIMD kernels
41  * which have explicit dependencies on the j-cluster size and/or SIMD-width.
42  * The functionality which depends on the j-cluster size is:
43  *   LJ-parameter lookup
44  *   force table lookup
45  *   energy group pair energy storage
46  */
47
48 /* Transpose 2 double precision registers */
49 static gmx_inline void
50 gmx_mm_transpose2_op_pd(__m128d in0, __m128d in1,
51                         __m128d *out0, __m128d *out1)
52 {
53     *out0 = _mm_unpacklo_pd(in0, in1);
54     *out1 = _mm_unpackhi_pd(in0, in1);
55 }
56
57 /* Sum the elements within each input register and store the sums in out */
58 static gmx_inline __m128d
59 gmx_mm_transpose_sum2_pr(__m128d in0, __m128d in1)
60 {
61     __m128d tr0, tr1;
62
63     gmx_mm_transpose2_op_pd(in0, in1, &tr0, &tr1);
64
65     return _mm_add_pd(tr0, tr1);
66 }
67
68 static inline __m128
69 gmx_mm128_invsqrt_ps_single(__m128 x)
70 {
71     const __m128 half  = _mm_set_ps(0.5, 0.5, 0.5, 0.5);
72     const __m128 three = _mm_set_ps(3.0, 3.0, 3.0, 3.0);
73
74     __m128       lu = _mm_rsqrt_ps(x);
75
76     return _mm_mul_ps(half, _mm_mul_ps(_mm_sub_ps(three, _mm_mul_ps(_mm_mul_ps(lu, lu), x)), lu));
77 }
78
79 /* Do 2 double precision invsqrt operations.
80  * Doing the SIMD rsqrt and the first Newton Raphson iteration
81  * in single precision gives full double precision accuracy.
82  */
83 static gmx_inline void
84 gmx_mm_invsqrt2_pd(__m128d in0, __m128d in1,
85                    __m128d *out0, __m128d *out1)
86 {
87     const __m128d half  = _mm_set1_pd(0.5);
88     const __m128d three = _mm_set1_pd(3.0);
89     __m128        s, ir;
90     __m128d       lu0, lu1;
91
92     s     = _mm_movelh_ps(_mm_cvtpd_ps(in0), _mm_cvtpd_ps(in1));
93     ir    = gmx_mm128_invsqrt_ps_single(s);
94     lu0   = _mm_cvtps_pd(ir);
95     lu1   = _mm_cvtps_pd(_mm_movehl_ps(ir, ir));
96     *out0 = _mm_mul_pd(half, _mm_mul_pd(_mm_sub_pd(three, _mm_mul_pd(_mm_mul_pd(lu0, lu0), in0)), lu0));
97     *out1 = _mm_mul_pd(half, _mm_mul_pd(_mm_sub_pd(three, _mm_mul_pd(_mm_mul_pd(lu1, lu1), in1)), lu1));
98 }
99
100 static gmx_inline void
101 load_lj_pair_params(const real *nbfp, const int *type, int aj,
102                     __m128d *c6_S, __m128d *c12_S)
103 {
104     __m128d clj_S[UNROLLJ];
105     int     p;
106
107     for (p = 0; p < UNROLLJ; p++)
108     {
109         clj_S[p] = _mm_load_pd(nbfp+type[aj+p]*NBFP_STRIDE);
110     }
111     gmx_mm_transpose2_op_pd(clj_S[0], clj_S[1], c6_S, c12_S);
112 }
113
114 /* The load_table functions below are performance critical.
115  * The routines issue UNROLLI*UNROLLJ _mm_load_ps calls.
116  * As these all have latencies, scheduling is crucial.
117  * The Intel compilers and CPUs seem to do a good job at this.
118  * But AMD CPUs perform significantly worse with gcc than with icc.
119  * Performance is improved a bit by using the extract function UNROLLJ times,
120  * instead of doing an _mm_store_si128 for every i-particle.
121  * This is only faster when we use FDV0 formatted tables, where we also need
122  * to multiple the index by 4, which can be done by a SIMD bit shift.
123  * With single precision AVX, 8 extracts are much slower than 1 store.
124  * Because of this, the load_table_f macro always takes the ti parameter,
125  * but it is only used with AVX.
126  */
127
128 static gmx_inline void
129 load_table_f(const real *tab_coul_F, gmx_epi32 ti_S, int *ti,
130              __m128d *ctab0_S, __m128d *ctab1_S)
131 {
132     int     idx[2];
133     __m128d ctab_S[2];
134
135     /* Without SSE4.1 the extract macro needs an immediate: unroll */
136     idx[0]    = gmx_mm_extract_epi32(ti_S, 0);
137     ctab_S[0] = _mm_loadu_pd(tab_coul_F+idx[0]);
138     idx[1]    = gmx_mm_extract_epi32(ti_S, 1);
139     ctab_S[1] = _mm_loadu_pd(tab_coul_F+idx[1]);
140
141     /* Shuffle the force table entries to a convenient order */
142     gmx_mm_transpose2_op_pd(ctab_S[0], ctab_S[1], ctab0_S, ctab1_S);
143     /* The second force table entry should contain the difference */
144     *ctab1_S = _mm_sub_pd(*ctab1_S, *ctab0_S);
145 }
146
147 static gmx_inline void
148 load_table_f_v(const real *tab_coul_F, const real *tab_coul_V,
149                gmx_epi32 ti_S, int *ti,
150                __m128d *ctab0_S, __m128d *ctab1_S, __m128d *ctabv_S)
151 {
152     int     idx[2];
153     __m128d ctab_S[4];
154
155     /* Without SSE4.1 the extract macro needs an immediate: unroll */
156     idx[0]    = gmx_mm_extract_epi32(ti_S, 0);
157     ctab_S[0] = _mm_loadu_pd(tab_coul_F+idx[0]);
158     idx[1]    = gmx_mm_extract_epi32(ti_S, 1);
159     ctab_S[1] = _mm_loadu_pd(tab_coul_F+idx[1]);
160
161     /* Shuffle the force table entries to a convenient order */
162     gmx_mm_transpose2_op_pd(ctab_S[0], ctab_S[1], ctab0_S, ctab1_S);
163     /* The second force table entry should contain the difference */
164     *ctab1_S = _mm_sub_pd(*ctab1_S, *ctab0_S);
165
166     ctab_S[2] = _mm_loadu_pd(tab_coul_V+idx[0]);
167     ctab_S[3] = _mm_loadu_pd(tab_coul_V+idx[1]);
168
169     /* Shuffle the energy table entries to a single register */
170     *ctabv_S = _mm_shuffle_pd(ctab_S[2], ctab_S[3], _MM_SHUFFLE2(0, 0));
171 }
172
173 #endif /* _nbnxn_kernel_simd_utils_x86_s128d_h_ */