SYCL: Use acc.bind(cgh) instead of cgh.require(acc)
[alexxy/gromacs.git] / src / gromacs / nbnxm / pairlistset.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014,2015,2017 by the GROMACS development team.
5  * Copyright (c) 2018,2019,2020,2021, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36
37 /*! \internal \file
38  *
39  * \brief
40  * Declares the PairlistSet class
41  *
42  * There is one PairlistSet object per locality. A PairlistSet
43  * holds a list of CPU- or GPU-type pairlist objects, one for each thread,
44  * as well as helper objects to construct each of those pairlists.
45  *
46  * \author Berk Hess <hess@kth.se>
47  * \ingroup module_nbnxm
48  */
49
50 #ifndef GMX_NBNXM_PAIRLISTSET_H
51 #define GMX_NBNXM_PAIRLISTSET_H
52
53 #include <memory>
54
55 #include "gromacs/math/vectypes.h"
56 #include "gromacs/mdtypes/locality.h"
57 #include "gromacs/utility/basedefinitions.h"
58 #include "gromacs/utility/real.h"
59
60 #include "pairlist.h"
61
62 struct nbnxn_atomdata_t;
63 struct PairlistParams;
64 struct PairsearchWork;
65 struct SearchCycleCounting;
66 struct t_nrnb;
67
68 namespace gmx
69 {
70 template<typename>
71 class ListOfLists;
72 }
73
74 namespace Nbnxm
75 {
76 class GridSet;
77 }
78
79 /*! \internal
80  * \brief An object that holds the local or non-local pairlists
81  */
82 class PairlistSet
83 {
84 public:
85     //! Constructor: initializes the pairlist set as empty
86     PairlistSet(const PairlistParams& listParams);
87
88     ~PairlistSet();
89
90     //! Constructs the pairlists in the set using the coordinates in \p nbat
91     void constructPairlists(gmx::InteractionLocality      locality,
92                             const Nbnxm::GridSet&         gridSet,
93                             gmx::ArrayRef<PairsearchWork> searchWork,
94                             nbnxn_atomdata_t*             nbat,
95                             const gmx::ListOfLists<int>&  exclusions,
96                             int                           minimumIlistCountForGpuBalancing,
97                             t_nrnb*                       nrnb,
98                             SearchCycleCounting*          searchCycleCounting);
99
100     //! Dispatch the kernel for dynamic pairlist pruning
101     void dispatchPruneKernel(const nbnxn_atomdata_t* nbat, gmx::ArrayRef<const gmx::RVec> shift_vec);
102
103     //! Returns the lists of CPU pairlists
104     gmx::ArrayRef<const NbnxnPairlistCpu> cpuLists() const { return cpuLists_; }
105
106     //! Returns a pointer to the GPU pairlist, nullptr when not present
107     const NbnxnPairlistGpu* gpuList() const
108     {
109         if (!gpuLists_.empty())
110         {
111             return &gpuLists_[0];
112         }
113         else
114         {
115             return nullptr;
116         }
117     }
118
119     //! Returns the lists of free-energy pairlists, empty when nonbonded interactions are not perturbed
120     gmx::ArrayRef<const std::unique_ptr<t_nblist>> fepLists() const { return fepLists_; }
121
122 private:
123     //! List of pairlists in CPU layout
124     std::vector<NbnxnPairlistCpu> cpuLists_;
125     //! List of working list for rebalancing CPU lists
126     std::vector<NbnxnPairlistCpu> cpuListsWork_;
127     //! List of pairlists in GPU layout
128     std::vector<NbnxnPairlistGpu> gpuLists_;
129     //! Pairlist parameters describing setup and ranges
130     const PairlistParams& params_;
131     //! Tells whether multiple lists get merged into one (the first) after creation
132     bool combineLists_;
133     //! Tells whether the lists is of CPU type, otherwise GPU type
134     gmx_bool isCpuType_;
135     //! Lists for perturbed interactions in simple atom-atom layout
136     std::vector<std::unique_ptr<t_nblist>> fepLists_;
137
138 public:
139     /* Pair counts for flop counting */
140     //! Total number of atom pairs for LJ+Q kernel
141     int natpair_ljq_;
142     //! Total number of atom pairs for LJ kernel
143     int natpair_lj_;
144     //! Total number of atom pairs for Q kernel
145     int natpair_q_;
146 };
147
148 #endif