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 the PairSearch class and helper structs
40 * The PairSearch class holds the domain setup, the search grids
41 * and helper object for the pair search. It manages the search work.
42 * The actual gridding and pairlist generation is performeed by the
43 * GridSet/Grid and PairlistSet/Pairlist classes, respectively.
45 * \author Berk Hess <hess@kth.se>
47 * \ingroup module_nbnxm
50 #ifndef GMX_NBNXM_PAIRSEARCH_H
51 #define GMX_NBNXM_PAIRSEARCH_H
56 #include "gromacs/domdec/domdec.h"
57 #include "gromacs/math/vectypes.h"
58 #include "gromacs/nbnxm/atomdata.h"
59 #include "gromacs/nbnxm/pairlist.h"
60 #include "gromacs/timing/cyclecounter.h"
61 #include "gromacs/utility/alignedallocator.h"
62 #include "gromacs/utility/arrayref.h"
63 #include "gromacs/utility/real.h"
67 struct gmx_domdec_zones_t
;
68 struct PairsearchWork
;
71 /*! \brief Convenience declaration for an std::vector with aligned memory */
73 using AlignedVector
= std::vector
< T
, gmx::AlignedAllocator
< T
>>;
76 /* Local cycle count struct for profiling */
82 start_
= gmx_cycles_read();
87 cycles_
+= gmx_cycles_read() - start_
;
96 double averageMCycles() const
100 return static_cast<double>(cycles_
)*1e-6/count_
;
110 gmx_cycles_t cycles_
= 0;
111 gmx_cycles_t start_
= 0;
114 //! Local cycle count enum for profiling different parts of search
116 enbsCCgrid
, enbsCCsearch
, enbsCCcombine
, enbsCCnr
120 * \brief Struct for collecting detailed cycle counts for the search
122 struct SearchCycleCounting
124 //! Start a pair search cycle counter
125 void start(const int enbsCC
)
130 //! Stop a pair search cycle counter
131 void stop(const int enbsCC
)
136 //! Print the cycle counts to \p fp
137 void printCycles(FILE *fp
,
138 gmx::ArrayRef
<const PairsearchWork
> work
) const;
140 //! Tells whether we record cycles
141 bool recordCycles_
= false;
142 //! The number of times pairsearching has been performed, local+non-local count as 1
143 int searchCount_
= 0;
144 //! The set of cycle counters
145 nbnxn_cycle_t cc_
[enbsCCnr
];
148 // TODO: Move nbnxn_search_work_t definition to its own file
150 /* Thread-local work struct, contains working data for Grid */
151 struct PairsearchWork
157 gmx_cache_protect_t cp0
; /* Buffer to avoid cache polution */
159 std::vector
<int> sortBuffer
; /* Temporary buffer for sorting atoms within a grid column */
161 nbnxn_buffer_flags_t buffer_flags
; /* Flags for force buffer access */
163 int ndistc
; /* Number of distance checks for flop counting */
166 std::unique_ptr
<t_nblist
> nbl_fep
; /* Temporary FEP list for load balancing */
168 nbnxn_cycle_t cycleCounter
; /* Counter for thread-local cycles */
170 gmx_cache_protect_t cp1
; /* Buffer to avoid cache polution */
173 /* Main pair-search struct, contains the grid(s), not the pair-list(s) */
177 //! Puts the atoms in \p ddZone on the grid and copies the coordinates to \p nbat
178 void putOnGrid(const matrix box
,
180 const rvec lowerCorner
,
181 const rvec upperCorner
,
182 const gmx::UpdateGroupsCog
*updateGroupsCog
,
186 gmx::ArrayRef
<const int> atomInfo
,
187 gmx::ArrayRef
<const gmx::RVec
> x
,
190 nbnxn_atomdata_t
*nbat
)
192 cycleCounting_
.start(enbsCCgrid
);
194 gridSet_
.putOnGrid(box
, ddZone
, lowerCorner
, upperCorner
,
195 updateGroupsCog
, atomStart
, atomEnd
, atomDensity
,
196 atomInfo
, x
, numAtomsMoved
, move
, nbat
);
198 cycleCounting_
.stop(enbsCCgrid
);
201 /* \brief Constructor
203 * \param[in] ePBC The periodic boundary conditions
204 * \param[in] numDDCells The number of domain decomposition cells per dimension, without DD nullptr should be passed
205 * \param[in] zones The domain decomposition zone setup, without DD nullptr should be passed
206 * \param[in] haveFep Tells whether non-bonded interactions are perturbed
207 * \param[in] maxNumThreads The maximum number of threads used in the search
210 const ivec
*numDDCells
,
211 const gmx_domdec_zones_t
*zones
,
212 PairlistType pairlistType
,
215 gmx::PinningPolicy pinningPolicy
);
217 //! Sets the order of the local atoms to the order grid atom ordering
218 void setLocalAtomOrder()
220 gridSet_
.setLocalAtomOrder();
223 //! Returns the set of search grids
224 const Nbnxm::GridSet
&gridSet() const
229 //! Returns the list of thread-local work objects
230 gmx::ArrayRef
<const PairsearchWork
> work() const
235 //! Returns the list of thread-local work objects
236 gmx::ArrayRef
<PairsearchWork
> work()
242 //! The set of search grids
243 Nbnxm::GridSet gridSet_
;
244 //! Work objects, one entry for each thread
245 std::vector
<PairsearchWork
> work_
;
248 //! Cycle counting for measuring components of the search
249 SearchCycleCounting cycleCounting_
;