9216457430e57b278fdd7818fc5202fa2cd03bda
[alexxy/gromacs.git] / src / gromacs / mdlib / nbnxn_kernels / simd_4xn / nbnxn_kernel_simd_4xn_common.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 #include "gromacs/pbcutil/ishift.h"
36 #include "gromacs/simd/simd.h"
37 #include "gromacs/simd/simd_math.h"
38 #include "gromacs/simd/vector_operations.h"
39 #include "../../nbnxn_consts.h"
40 #ifdef CALC_COUL_EWALD
41 #include "gromacs/math/utilities.h"
42 #endif
43
44 #include "config.h"
45
46 #ifndef GMX_SIMD_J_UNROLL_SIZE
47 #error "Need to define GMX_SIMD_J_UNROLL_SIZE before including the 4xn kernel common header file"
48 #endif
49
50 #define UNROLLI    NBNXN_CPU_CLUSTER_I_SIZE
51 #define UNROLLJ    (GMX_SIMD_REAL_WIDTH/GMX_SIMD_J_UNROLL_SIZE)
52
53 /* The stride of all the atom data arrays is max(UNROLLI,unrollj) */
54 #if GMX_SIMD_REAL_WIDTH >= UNROLLI
55 #define STRIDE     (GMX_SIMD_REAL_WIDTH/GMX_SIMD_J_UNROLL_SIZE)
56 #else
57 #define STRIDE     (UNROLLI)
58 #endif
59
60 #include "../nbnxn_kernel_simd_utils.h"
61
62 static gmx_inline void gmx_simdcall
63 gmx_load_simd_4xn_interactions(int gmx_unused             excl,
64                                gmx_exclfilter gmx_unused  filter_S0,
65                                gmx_exclfilter gmx_unused  filter_S1,
66                                gmx_exclfilter gmx_unused  filter_S2,
67                                gmx_exclfilter gmx_unused  filter_S3,
68                                const char gmx_unused     *interaction_mask_indices,
69                                real gmx_unused           *simd_interaction_array,
70                                gmx_simd_bool_t           *interact_S0,
71                                gmx_simd_bool_t           *interact_S1,
72                                gmx_simd_bool_t           *interact_S2,
73                                gmx_simd_bool_t           *interact_S3)
74 {
75 #if defined GMX_SIMD_X86_SSE2_OR_HIGHER || defined GMX_SIMD_REFERENCE
76     /* Load integer interaction mask */
77     gmx_exclfilter mask_pr_S = gmx_load1_exclfilter(excl);
78     *interact_S0  = gmx_checkbitmask_pb(mask_pr_S, filter_S0);
79     *interact_S1  = gmx_checkbitmask_pb(mask_pr_S, filter_S1);
80     *interact_S2  = gmx_checkbitmask_pb(mask_pr_S, filter_S2);
81     *interact_S3  = gmx_checkbitmask_pb(mask_pr_S, filter_S3);
82 #endif
83 #ifdef GMX_SIMD_IBM_QPX
84     const int size = GMX_SIMD_REAL_WIDTH * sizeof(real);
85     *interact_S0  = gmx_load_interaction_mask_pb(size*interaction_mask_indices[0], simd_interaction_array);
86     *interact_S1  = gmx_load_interaction_mask_pb(size*interaction_mask_indices[1], simd_interaction_array);
87     *interact_S2  = gmx_load_interaction_mask_pb(size*interaction_mask_indices[2], simd_interaction_array);
88     *interact_S3  = gmx_load_interaction_mask_pb(size*interaction_mask_indices[3], simd_interaction_array);
89 #endif
90 }
91
92 /* All functionality defines are set here, except for:
93  * CALC_ENERGIES, ENERGY_GROUPS which are defined before.
94  * CHECK_EXCLS, which is set just before including the inner loop contents.
95  * The combination rule defines, LJ_COMB_GEOM or LJ_COMB_LB are currently
96  * set before calling the kernel function. We might want to move that
97  * to inside the n-loop and have a different combination rule for different
98  * ci's, as no combination rule gives a 50% performance hit for LJ.
99  */
100
101 /* We always calculate shift forces, because it's cheap anyhow */
102 #define CALC_SHIFTFORCES
103
104 /* Assumes all LJ parameters are identical */
105 /* #define FIX_LJ_C */