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/atomdata.h"
54 #include "gromacs/nbnxm/pairlist.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"
62 struct gmx_domdec_zones_t;
65 /*! \brief Convenience declaration for an std::vector with aligned memory */
67 using AlignedVector = std::vector < T, gmx::AlignedAllocator < T>>;
70 /* Local cycle count struct for profiling */
76 start_ = gmx_cycles_read();
81 cycles_ += gmx_cycles_read() - start_;
90 double averageMCycles() const
94 return static_cast<double>(cycles_)*1e-6/count_;
104 gmx_cycles_t cycles_ = 0;
105 gmx_cycles_t start_ = 0;
108 // TODO: Move nbnxn_search_work_t definition to its own file
110 /* Thread-local work struct, contains working data for Grid */
111 struct PairsearchWork
117 gmx_cache_protect_t cp0; /* Buffer to avoid cache polution */
119 std::vector<int> sortBuffer; /* Temporary buffer for sorting atoms within a grid column */
121 nbnxn_buffer_flags_t buffer_flags; /* Flags for force buffer access */
123 int ndistc; /* Number of distance checks for flop counting */
126 std::unique_ptr<t_nblist> nbl_fep; /* Temporary FEP list for load balancing */
128 nbnxn_cycle_t cycleCounter; /* Counter for thread-local cycles */
130 gmx_cache_protect_t cp1; /* Buffer to avoid cache polution */
133 /* Main pair-search struct, contains the grid(s), not the pair-list(s) */
138 * \brief Description of the domain setup: PBC and the connections between domains
143 * \brief Description of the domain setup: PBC and the connections between domains
145 //! Constructor, without DD \p numDDCells and \p ddZones should be nullptr
146 DomainSetup(int ePBC,
147 const ivec *numDDCells,
148 const gmx_domdec_zones_t *ddZones);
152 //! Tells whether we are using domain decomposition
154 //! Tells whether we are using domain decomposition per dimension
155 std::array<bool, DIM> haveDomDecPerDim;
156 //! The domain decomposition zone setup
157 const gmx_domdec_zones_t *zones;
160 //! Local cycle count enum for profiling different parts of search
162 enbsCCgrid, enbsCCsearch, enbsCCcombine, enbsCCnr
165 struct SearchCycleCounting
167 //! Start a pair search cycle counter
168 void start(const int enbsCC)
173 //! Stop a pair search cycle counter
174 void stop(const int enbsCC)
179 //! Print the cycle counts to \p fp
180 void printCycles(FILE *fp,
181 gmx::ArrayRef<const PairsearchWork> work) const;
183 bool recordCycles_ = false;
184 int searchCount_ = 0;
185 nbnxn_cycle_t cc_[enbsCCnr];
188 //! Puts the atoms in \p ddZone on the grid and copies the coordinates to \p nbat
189 void putOnGrid(const matrix box,
191 const rvec lowerCorner,
192 const rvec upperCorner,
193 const gmx::UpdateGroupsCog *updateGroupsCog,
198 gmx::ArrayRef<const gmx::RVec> x,
201 nbnxn_atomdata_t *nbat)
203 cycleCounting_.start(enbsCCgrid);
205 gridSet_.putOnGrid(box, ddZone, lowerCorner, upperCorner,
206 updateGroupsCog, atomStart, atomEnd, atomDensity,
207 atinfo, x, numAtomsMoved, move, nbat);
209 cycleCounting_.stop(enbsCCgrid);
212 /* \brief Constructor
214 * \param[in] ePBC The periodic boundary conditions
215 * \param[in] numDDCells The number of domain decomposition cells per dimension, without DD nullptr should be passed
216 * \param[in] zones The domain decomposition zone setup, without DD nullptr should be passed
217 * \param[in] haveFep Tells whether non-bonded interactions are perturbed
218 * \param[in] maxNumThreads The maximum number of threads used in the search
221 const ivec *numDDCells,
222 const gmx_domdec_zones_t *zones,
223 PairlistType pairlistType,
227 //! Sets the order of the local atoms to the order grid atom ordering
228 void setLocalAtomOrder()
230 gridSet_.setLocalAtomOrder();
233 const DomainSetup domainSetup() const
238 //! Returns the set of search grids
239 const Nbnxm::GridSet &gridSet() const
244 //! Returns the list of thread-local work objects
245 gmx::ArrayRef<const PairsearchWork> work() const
250 //! Returns the list of thread-local work objects
251 gmx::ArrayRef<PairsearchWork> work()
258 DomainSetup domainSetup_;
259 //! The set of search grids
260 Nbnxm::GridSet gridSet_;
261 //! Work objects, one entry for each thread
262 std::vector<PairsearchWork> work_;
265 //! Cycle counting for measuring components of the search
266 SearchCycleCounting cycleCounting_;