Move PairlistSet declaration
[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,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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35
36 /*! \internal \file
37  *
38  * \brief
39  * Declares the PairlistSet class
40  *
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.
44  *
45  * \author Berk Hess <hess@kth.se>
46  * \ingroup module_nbnxm
47  */
48
49 #ifndef GMX_NBNXM_PAIRLISTSET_H
50 #define GMX_NBNXM_PAIRLISTSET_H
51
52 #include "gromacs/math/vectypes.h"
53 #include "gromacs/nbnxm/nbnxm.h"
54 #include "gromacs/nbnxm/pairlist.h"
55 #include "gromacs/utility/basedefinitions.h"
56 #include "gromacs/utility/real.h"
57
58 #include "locality.h"
59
60 struct PairlistParams;
61 struct PairsearchWork;
62 struct SearchCycleCounting;
63 struct t_nrnb;
64
65 namespace Nbnxm
66 {
67 class GridSet;
68 }
69
70 /*! \internal
71  * \brief An object that holds the local or non-local pairlists
72  */
73 class PairlistSet
74 {
75     public:
76         //! Constructor: initializes the pairlist set as empty
77         PairlistSet(Nbnxm::InteractionLocality  locality,
78                     const PairlistParams       &listParams);
79
80         ~PairlistSet();
81
82         //! Constructs the pairlists in the set using the coordinates in \p nbat
83         void constructPairlists(const Nbnxm::GridSet          &gridSet,
84                                 gmx::ArrayRef<PairsearchWork>  searchWork,
85                                 nbnxn_atomdata_t              *nbat,
86                                 const t_blocka                *excl,
87                                 Nbnxm::KernelType              kernelType,
88                                 int                            minimumIlistCountForGpuBalancing,
89                                 t_nrnb                        *nrnb,
90                                 SearchCycleCounting           *searchCycleCounting);
91
92         //! Dispatch the kernel for dynamic pairlist pruning
93         void dispatchPruneKernel(const nbnxn_atomdata_t *nbat,
94                                  const rvec             *shift_vec,
95                                  Nbnxm::KernelType       kernelType);
96
97         //! Returns the locality
98         Nbnxm::InteractionLocality locality() const
99         {
100             return locality_;
101         }
102
103         //! Returns the lists of CPU pairlists
104         gmx::ArrayRef<const NbnxnPairlistCpu> cpuLists() const
105         {
106             return cpuLists_;
107         }
108
109         //! Returns a pointer to the GPU pairlist, nullptr when not present
110         const NbnxnPairlistGpu *gpuList() const
111         {
112             if (!gpuLists_.empty())
113             {
114                 return &gpuLists_[0];
115             }
116             else
117             {
118                 return nullptr;
119             }
120         }
121
122         //! Returns the lists of free-energy pairlists, empty when nonbonded interactions are not perturbed
123         gmx::ArrayRef<t_nblist const * const> fepLists() const
124         {
125             return fepLists_;
126         }
127
128     private:
129         //! The locality of the pairlist set
130         Nbnxm::InteractionLocality     locality_;
131         //! List of pairlists in CPU layout
132         std::vector<NbnxnPairlistCpu>  cpuLists_;
133         //! List of working list for rebalancing CPU lists
134         std::vector<NbnxnPairlistCpu>  cpuListsWork_;
135         //! List of pairlists in GPU layout
136         std::vector<NbnxnPairlistGpu>  gpuLists_;
137         //! Pairlist parameters describing setup and ranges
138         const PairlistParams          &params_;
139         //! Tells whether multiple lists get merged into one (the first) after creation
140         bool                           combineLists_;
141         //! Tells whether the lists is of CPU type, otherwise GPU type
142         gmx_bool                       isCpuType_;
143         //! Lists for perturbed interactions in simple atom-atom layout
144         std::vector<t_nblist *>        fepLists_;
145
146     public:
147         /* Pair counts for flop counting */
148         //! Total number of atom pairs for LJ+Q kernel
149         int natpair_ljq_;
150         //! Total number of atom pairs for LJ kernel
151         int natpair_lj_;
152         //! Total number of atom pairs for Q kernel
153         int natpair_q_;
154 };
155
156 #endif