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/nbnxm/pairlist.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
;
74 * \brief An object that holds the local or non-local pairlists
79 //! Constructor: initializes the pairlist set as empty
80 PairlistSet(Nbnxm::InteractionLocality locality
,
81 const PairlistParams
&listParams
);
85 //! Constructs the pairlists in the set using the coordinates in \p nbat
86 void constructPairlists(const Nbnxm::GridSet
&gridSet
,
87 gmx::ArrayRef
<PairsearchWork
> searchWork
,
88 nbnxn_atomdata_t
*nbat
,
90 int minimumIlistCountForGpuBalancing
,
92 SearchCycleCounting
*searchCycleCounting
);
94 //! Dispatch the kernel for dynamic pairlist pruning
95 void dispatchPruneKernel(const nbnxn_atomdata_t
*nbat
,
96 const rvec
*shift_vec
);
98 //! Returns the locality
99 Nbnxm::InteractionLocality
locality() const
104 //! Returns the lists of CPU pairlists
105 gmx::ArrayRef
<const NbnxnPairlistCpu
> cpuLists() const
110 //! Returns a pointer to the GPU pairlist, nullptr when not present
111 const NbnxnPairlistGpu
*gpuList() const
113 if (!gpuLists_
.empty())
115 return &gpuLists_
[0];
123 //! Returns the lists of free-energy pairlists, empty when nonbonded interactions are not perturbed
124 gmx::ArrayRef
< const std::unique_ptr
< t_nblist
>> fepLists() const
130 //! The locality of the pairlist set
131 Nbnxm::InteractionLocality locality_
;
132 //! List of pairlists in CPU layout
133 std::vector
<NbnxnPairlistCpu
> cpuLists_
;
134 //! List of working list for rebalancing CPU lists
135 std::vector
<NbnxnPairlistCpu
> cpuListsWork_
;
136 //! List of pairlists in GPU layout
137 std::vector
<NbnxnPairlistGpu
> gpuLists_
;
138 //! Pairlist parameters describing setup and ranges
139 const PairlistParams
¶ms_
;
140 //! Tells whether multiple lists get merged into one (the first) after creation
142 //! Tells whether the lists is of CPU type, otherwise GPU type
144 //! Lists for perturbed interactions in simple atom-atom layout
145 std::vector
< std::unique_ptr
< t_nblist
>> fepLists_
;
148 /* Pair counts for flop counting */
149 //! Total number of atom pairs for LJ+Q kernel
151 //! Total number of atom pairs for LJ kernel
153 //! Total number of atom pairs for Q kernel