1 /*
2 * This file is part of the GROMACS molecular simulation package.
3 *
4 * Copyright (c) 2012,2013,2014,2015,2018 by the GROMACS development team.
5 * Copyright (c) 2019,2020, by the GROMACS development team, led by
6 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7 * and including many others, as listed in the AUTHORS file in the
8 * top-level source directory and at http://www.gromacs.org.
9 *
10 * GROMACS is free software; you can redistribute it and/or
11 * modify it under the terms of the GNU Lesser General Public License
12 * as published by the Free Software Foundation; either version 2.1
13 * of the License, or (at your option) any later version.
14 *
15 * GROMACS is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18 * Lesser General Public License for more details.
19 *
20 * You should have received a copy of the GNU Lesser General Public
21 * License along with GROMACS; if not, see
22 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 *
25 * If you want to redistribute modifications to GROMACS, please
26 * consider that scientific software is very special. Version
27 * control is crucial - bugs must be traceable. We will be happy to
28 * consider code for inclusion in the official distribution, but
29 * derived work must not be called official GROMACS. Details are found
30 * in the README & COPYING files - if they are missing, get the
31 * official version at http://www.gromacs.org.
32 *
33 * To help us fund GROMACS development, we humbly ask that you cite
34 * the research papers on the package. Check out http://www.gromacs.org.
35 */
36 #include "gromacs/pbcutil/ishift.h"
37 #include "gromacs/simd/simd.h"
38 #include "gromacs/simd/simd_math.h"
39 #include "gromacs/simd/vector_operations.h"
40 #ifdef CALC_COUL_EWALD
41 # include "gromacs/math/utilities.h"
42 #endif
43
44 #include "config.h"
45
46 #include <cstdint>
47
48 #if !GMX_SIMD_HAVE_HSIMD_UTIL_REAL
49 # error "Half-simd utility operations are required for the 2xNN kernels"
50 #endif
51
52 #ifndef GMX_SIMD_J_UNROLL_SIZE
53 # error "Need to define GMX_SIMD_J_UNROLL_SIZE before including the 2xnn kernel common header file"
54 #endif
55
56 #define UNROLLI 4
57 #define UNROLLJ (GMX_SIMD_REAL_WIDTH / GMX_SIMD_J_UNROLL_SIZE)
58
59 static_assert(UNROLLI == c_nbnxnCpuIClusterSize, "UNROLLI should match the i-cluster size");
60
61 /* The stride of all the atom data arrays is equal to half the SIMD width */
62 #define STRIDE UNROLLJ
63
64 #if !defined GMX_NBNXN_SIMD_2XNN && !defined GMX_NBNXN_SIMD_4XN
65 # error "Must define an NBNxN kernel flavour before including NBNxN kernel utility functions"
66 #endif
67
68 // We use the FDV0 tables for width==4 (when we can load it in one go), or if we don't have any unaligned loads
69 #if GMX_SIMD_REAL_WIDTH == 4 || !GMX_SIMD_HAVE_GATHER_LOADU_BYSIMDINT_TRANSPOSE_REAL
70 # define TAB_FDV0
71 #endif
72
73 #if defined UNROLLJ
74 /* As add_ener_grp, but for two groups of UNROLLJ/2 stored in
75 * a single SIMD register.
76 */
add_ener_grp_halves(gmx::SimdReal e_S,real * v0,real * v1,const int * offset_jj)77 static inline void add_ener_grp_halves(gmx::SimdReal e_S, real* v0, real* v1, const int* offset_jj)
78 {
79 for (int jj = 0; jj < (UNROLLJ / 2); jj++)
80 {
81 incrDualHsimd(v0 + offset_jj[jj] + jj * GMX_SIMD_REAL_WIDTH / 2,
82 v1 + offset_jj[jj] + jj * GMX_SIMD_REAL_WIDTH / 2, e_S);
83 }
84 }
85 #endif
86
87 #if GMX_SIMD_HAVE_INT32_LOGICAL
88 typedef gmx::SimdInt32 SimdBitMask;
89 #else
90 typedef gmx::SimdReal SimdBitMask;
91 #endif
92
93
gmx_load_simd_2xnn_interactions(int excl,SimdBitMask filter_S0,SimdBitMask filter_S2,gmx::SimdBool * interact_S0,gmx::SimdBool * interact_S2)94 static inline void gmx_simdcall gmx_load_simd_2xnn_interactions(int excl,
95 SimdBitMask filter_S0,
96 SimdBitMask filter_S2,
97 gmx::SimdBool* interact_S0,
98 gmx::SimdBool* interact_S2)
99 {
100 using namespace gmx;
101 #if GMX_SIMD_HAVE_INT32_LOGICAL
102 SimdInt32 mask_pr_S(excl);
103 *interact_S0 = cvtIB2B(testBits(mask_pr_S & filter_S0));
104 *interact_S2 = cvtIB2B(testBits(mask_pr_S & filter_S2));
105 #elif GMX_SIMD_HAVE_LOGICAL
106 union {
107 # if GMX_DOUBLE
108 std::int64_t i;
109 # else
110 std::int32_t i;
111 # endif
112 real r;
113 } conv;
114
115 conv.i = excl;
116 SimdReal mask_pr_S(conv.r);
117
118 *interact_S0 = testBits(mask_pr_S & filter_S0);
119 *interact_S2 = testBits(mask_pr_S & filter_S2);
120 #endif
121 }
122
123 /* All functionality defines are set here, except for:
124 * CALC_ENERGIES, ENERGY_GROUPS which are defined before.
125 * CHECK_EXCLS, which is set just before including the inner loop contents.
126 * The combination rule defines, LJ_COMB_GEOM or LJ_COMB_LB are currently
127 * set before calling the kernel function. We might want to move that
128 * to inside the n-loop and have a different combination rule for different
129 * ci's, as no combination rule gives a 50% performance hit for LJ.
130 */
131
132 /* We always calculate shift forces, because it's cheap anyhow */
133 #define CALC_SHIFTFORCES
134
135 /* Assumes all LJ parameters are identical */
136 /* #define FIX_LJ_C */
137