2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,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.
39 * Declares the PairlistSet class
41 * There is one PairlistSet object per locality. A PairlistSet
42 * holds a list of CPU- or GPU-type pairlist objects, one for each thread,
43 * as well as helper objects to construct each of those pairlists.
45 * \author Berk Hess <hess@kth.se>
46 * \ingroup module_nbnxm
49 #ifndef GMX_NBNXM_PAIRLISTSET_H
50 #define GMX_NBNXM_PAIRLISTSET_H
54 #include "gromacs/math/vectypes.h"
55 #include "gromacs/mdtypes/locality.h"
56 #include "gromacs/utility/basedefinitions.h"
57 #include "gromacs/utility/real.h"
61 struct nbnxn_atomdata_t;
62 struct PairlistParams;
63 struct PairsearchWork;
64 struct SearchCycleCounting;
79 * \brief An object that holds the local or non-local pairlists
84 //! Constructor: initializes the pairlist set as empty
85 PairlistSet(gmx::InteractionLocality locality, const PairlistParams& listParams);
89 //! Constructs the pairlists in the set using the coordinates in \p nbat
90 void constructPairlists(const Nbnxm::GridSet& gridSet,
91 gmx::ArrayRef<PairsearchWork> searchWork,
92 nbnxn_atomdata_t* nbat,
93 const gmx::ListOfLists<int>& exclusions,
94 int minimumIlistCountForGpuBalancing,
96 SearchCycleCounting* searchCycleCounting);
98 //! Dispatch the kernel for dynamic pairlist pruning
99 void dispatchPruneKernel(const nbnxn_atomdata_t* nbat, const rvec* shift_vec);
101 //! Returns the locality
102 gmx::InteractionLocality locality() const { return locality_; }
104 //! Returns the lists of CPU pairlists
105 gmx::ArrayRef<const NbnxnPairlistCpu> cpuLists() const { return cpuLists_; }
107 //! Returns a pointer to the GPU pairlist, nullptr when not present
108 const NbnxnPairlistGpu* gpuList() const
110 if (!gpuLists_.empty())
112 return &gpuLists_[0];
120 //! Returns the lists of free-energy pairlists, empty when nonbonded interactions are not perturbed
121 gmx::ArrayRef<const std::unique_ptr<t_nblist>> fepLists() const { return fepLists_; }
124 //! The locality of the pairlist set
125 gmx::InteractionLocality locality_;
126 //! List of pairlists in CPU layout
127 std::vector<NbnxnPairlistCpu> cpuLists_;
128 //! List of working list for rebalancing CPU lists
129 std::vector<NbnxnPairlistCpu> cpuListsWork_;
130 //! List of pairlists in GPU layout
131 std::vector<NbnxnPairlistGpu> gpuLists_;
132 //! Pairlist parameters describing setup and ranges
133 const PairlistParams& params_;
134 //! Tells whether multiple lists get merged into one (the first) after creation
136 //! Tells whether the lists is of CPU type, otherwise GPU type
138 //! Lists for perturbed interactions in simple atom-atom layout
139 std::vector<std::unique_ptr<t_nblist>> fepLists_;
142 /* Pair counts for flop counting */
143 //! Total number of atom pairs for LJ+Q kernel
145 //! Total number of atom pairs for LJ kernel
147 //! Total number of atom pairs for Q kernel