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