2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016,2017,2018,2019, 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.
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.
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.
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.
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.
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.
38 * \brief Declares internal nbnxm module details
40 * \author Berk Hess <hess@kth.se>
42 * \ingroup module_nbnxm
45 #ifndef GMX_NBNXM_INTERNAL_H
46 #define GMX_NBNXM_INTERNAL_H
51 #include "gromacs/domdec/domdec.h"
52 #include "gromacs/math/vectypes.h"
53 #include "gromacs/nbnxm/pairlist.h"
54 #include "gromacs/simd/simd.h"
55 #include "gromacs/timing/cyclecounter.h"
56 #include "gromacs/utility/alignedallocator.h"
57 #include "gromacs/utility/arrayref.h"
58 #include "gromacs/utility/real.h"
60 struct gmx_domdec_zones_t;
68 // TODO Document after refactoring
71 /* Strides for x/f with xyz and xyzq coordinate (and charge) storage */
74 /* Size of packs of x, y or z with SIMD packed coords/forces */
75 static const int c_packX4 = 4;
76 static const int c_packX8 = 8;
77 /* Strides for a pack of 4 and 8 coordinates/forces */
78 #define STRIDE_P4 (DIM*c_packX4)
79 #define STRIDE_P8 (DIM*c_packX8)
81 /* Returns the index in a coordinate array corresponding to atom a */
82 template<int packSize> static inline int atom_to_x_index(int a)
84 return DIM*(a & ~(packSize - 1)) + (a & (packSize - 1));
89 /* Memory alignment in bytes as required by SIMD aligned loads/stores */
90 #define NBNXN_MEM_ALIGN (GMX_SIMD_REAL_WIDTH*sizeof(real))
92 /* No alignment required, but set it so we can call the same routines */
93 #define NBNXN_MEM_ALIGN 32
99 /*! \brief Convenience declaration for an std::vector with aligned memory */
101 using AlignedVector = std::vector < T, gmx::AlignedAllocator < T>>;
104 /* Local cycle count struct for profiling */
111 /* Local cycle count enum for profiling */
113 enbsCCgrid, enbsCCsearch, enbsCCcombine, enbsCCreducef, enbsCCnr
116 /* Thread-local work struct, contains working data for Grid */
117 struct nbnxn_search_work_t
119 nbnxn_search_work_t();
121 ~nbnxn_search_work_t();
123 gmx_cache_protect_t cp0; /* Buffer to avoid cache polution */
125 std::vector<int> cxy_na; /* Grid column atom counts temporary buffer */
127 std::vector<int> sortBuffer; /* Temporary buffer for sorting atoms within a grid column */
129 nbnxn_buffer_flags_t buffer_flags; /* Flags for force buffer access */
131 int ndistc; /* Number of distance checks for flop counting */
134 std::unique_ptr<t_nblist> nbl_fep; /* Temporary FEP list for load balancing */
136 nbnxn_cycle_t cc[enbsCCnr]; /* Counters for thread-local cycles */
138 gmx_cache_protect_t cp1; /* Buffer to avoid cache polution */
141 /* Main pair-search struct, contains the grid(s), not the pair-list(s) */
144 /* \brief Constructor
146 * \param[in] ePBC The periodic boundary conditions
147 * \param[in] n_dd_cells The number of domain decomposition cells per dimension, without DD nullptr should be passed
148 * \param[in] zones The domain decomposition zone setup, without DD nullptr should be passed
149 * \param[in] nb_kernel_type The nbnxn non-bonded kernel type
150 * \param[in] bFEP Tells whether non-bonded interactions are perturbed
151 * \param[in] nthread_max The maximum number of threads used in the search
154 nbnxn_search(int ePBC,
155 const ivec *n_dd_cells,
156 const gmx_domdec_zones_t *zones,
157 PairlistType pairlistType,
161 gmx_bool bFEP; /* Do we have perturbed atoms? */
162 int ePBC; /* PBC type enum */
163 matrix box; /* The periodic unit-cell */
165 gmx_bool DomDec; /* Are we doing domain decomposition? */
166 ivec dd_dim; /* Are we doing DD in x,y,z? */
167 const gmx_domdec_zones_t *zones; /* The domain decomposition zones */
169 std::vector<Nbnxm::Grid> grid; /* Array of grids, size ngrid */
170 std::vector<int> cell; /* Actual allocated cell array for all grids */
171 std::vector<int> a; /* Atom index for grid, the inverse of cell */
173 int natoms_local; /* The local atoms run from 0 to natoms_local */
174 int natoms_nonlocal; /* The non-local atoms run from natoms_local
175 * to natoms_nonlocal */
177 gmx_bool print_cycles;
179 nbnxn_cycle_t cc[enbsCCnr];
181 /* Thread-local work data */
182 mutable std::vector<nbnxn_search_work_t> work; /* Work array, one entry for each thread */
186 /*! \brief Start an nbnxn cycle counter */
187 static inline void nbs_cycle_start(nbnxn_cycle_t *cc)
189 cc->start = gmx_cycles_read();
192 /*! \brief Stop an nbnxn cycle counter */
193 static inline void nbs_cycle_stop(nbnxn_cycle_t *cc)
195 cc->c += gmx_cycles_read() - cc->start;