Apply clang-format to source tree
[alexxy/gromacs.git] / src / gromacs / nbnxm / pairsearch.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014,2015,2016,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 Declares the PairSearch class and helper structs
39  *
40  * The PairSearch class holds the domain setup, the search grids
41  * and helper object for the pair search. It manages the search work.
42  * The actual gridding and pairlist generation is performeed by the
43  * GridSet/Grid and PairlistSet/Pairlist classes, respectively.
44  *
45  * \author Berk Hess <hess@kth.se>
46  *
47  * \ingroup module_nbnxm
48  */
49
50 #ifndef GMX_NBNXM_PAIRSEARCH_H
51 #define GMX_NBNXM_PAIRSEARCH_H
52
53 #include <memory>
54 #include <vector>
55
56 #include "gromacs/domdec/domdec.h"
57 #include "gromacs/math/vectypes.h"
58 #include "gromacs/timing/cyclecounter.h"
59 #include "gromacs/utility/alignedallocator.h"
60 #include "gromacs/utility/arrayref.h"
61 #include "gromacs/utility/real.h"
62
63 #include "atomdata.h"
64 #include "gridset.h"
65 #include "pairlist.h"
66
67 struct gmx_domdec_zones_t;
68 struct PairsearchWork;
69
70
71 /*! \brief Convenience declaration for an std::vector with aligned memory */
72 template<class T>
73 using AlignedVector = std::vector<T, gmx::AlignedAllocator<T>>;
74
75
76 /* Local cycle count struct for profiling */
77 class nbnxn_cycle_t
78 {
79 public:
80     void start() { start_ = gmx_cycles_read(); }
81
82     void stop()
83     {
84         cycles_ += gmx_cycles_read() - start_;
85         count_++;
86     }
87
88     int count() const { return count_; }
89
90     double averageMCycles() const
91     {
92         if (count_ > 0)
93         {
94             return static_cast<double>(cycles_) * 1e-6 / count_;
95         }
96         else
97         {
98             return 0;
99         }
100     }
101
102 private:
103     int          count_  = 0;
104     gmx_cycles_t cycles_ = 0;
105     gmx_cycles_t start_  = 0;
106 };
107
108 //! Local cycle count enum for profiling different parts of search
109 enum
110 {
111     enbsCCgrid,
112     enbsCCsearch,
113     enbsCCcombine,
114     enbsCCnr
115 };
116
117 /*! \internal
118  * \brief Struct for collecting detailed cycle counts for the search
119  */
120 struct SearchCycleCounting
121 {
122     //! Start a pair search cycle counter
123     void start(const int enbsCC) { cc_[enbsCC].start(); }
124
125     //! Stop a pair search cycle counter
126     void stop(const int enbsCC) { cc_[enbsCC].stop(); }
127
128     //! Print the cycle counts to \p fp
129     void printCycles(FILE* fp, gmx::ArrayRef<const PairsearchWork> work) const;
130
131     //! Tells whether we record cycles
132     bool recordCycles_ = false;
133     //! The number of times pairsearching has been performed, local+non-local count as 1
134     int searchCount_ = 0;
135     //! The set of cycle counters
136     nbnxn_cycle_t cc_[enbsCCnr];
137 };
138
139 // TODO: Move nbnxn_search_work_t definition to its own file
140
141 /* Thread-local work struct, contains working data for Grid */
142 struct PairsearchWork
143 {
144     PairsearchWork();
145
146     ~PairsearchWork();
147
148     gmx_cache_protect_t cp0; /* Buffer to avoid cache polution */
149
150     std::vector<int> sortBuffer; /* Temporary buffer for sorting atoms within a grid column */
151
152     nbnxn_buffer_flags_t buffer_flags; /* Flags for force buffer access */
153
154     int ndistc; /* Number of distance checks for flop counting */
155
156
157     std::unique_ptr<t_nblist> nbl_fep; /* Temporary FEP list for load balancing */
158
159     nbnxn_cycle_t cycleCounter; /* Counter for thread-local cycles */
160
161     gmx_cache_protect_t cp1; /* Buffer to avoid cache polution */
162 };
163
164 /* Main pair-search struct, contains the grid(s), not the pair-list(s) */
165 class PairSearch
166 {
167 public:
168     //! Puts the atoms in \p ddZone on the grid and copies the coordinates to \p nbat
169     void putOnGrid(const matrix                   box,
170                    int                            ddZone,
171                    const rvec                     lowerCorner,
172                    const rvec                     upperCorner,
173                    const gmx::UpdateGroupsCog*    updateGroupsCog,
174                    gmx::Range<int>                atomRange,
175                    real                           atomDensity,
176                    gmx::ArrayRef<const int>       atomInfo,
177                    gmx::ArrayRef<const gmx::RVec> x,
178                    int                            numAtomsMoved,
179                    const int*                     move,
180                    nbnxn_atomdata_t*              nbat)
181     {
182         cycleCounting_.start(enbsCCgrid);
183
184         gridSet_.putOnGrid(box, ddZone, lowerCorner, upperCorner, updateGroupsCog, atomRange,
185                            atomDensity, atomInfo, x, numAtomsMoved, move, nbat);
186
187         cycleCounting_.stop(enbsCCgrid);
188     }
189
190     /* \brief Constructor
191      *
192      * \param[in] ePBC            The periodic boundary conditions
193      * \param[in] numDDCells      The number of domain decomposition cells per dimension, without DD nullptr should be passed
194      * \param[in] zones           The domain decomposition zone setup, without DD nullptr should be passed
195      * \param[in] haveFep         Tells whether non-bonded interactions are perturbed
196      * \param[in] maxNumThreads   The maximum number of threads used in the search
197      */
198     PairSearch(int                       ePBC,
199                bool                      doTestParticleInsertion,
200                const ivec*               numDDCells,
201                const gmx_domdec_zones_t* zones,
202                PairlistType              pairlistType,
203                bool                      haveFep,
204                int                       maxNumthreads,
205                gmx::PinningPolicy        pinningPolicy);
206
207     //! Sets the order of the local atoms to the order grid atom ordering
208     void setLocalAtomOrder() { gridSet_.setLocalAtomOrder(); }
209
210     //! Returns the set of search grids
211     const Nbnxm::GridSet& gridSet() const { return gridSet_; }
212
213     //! Returns the list of thread-local work objects
214     gmx::ArrayRef<const PairsearchWork> work() const { return work_; }
215
216     //! Returns the list of thread-local work objects
217     gmx::ArrayRef<PairsearchWork> work() { return work_; }
218
219 private:
220     //! The set of search grids
221     Nbnxm::GridSet gridSet_;
222     //! Work objects, one entry for each thread
223     std::vector<PairsearchWork> work_;
224
225 public:
226     //! Cycle counting for measuring components of the search
227     SearchCycleCounting cycleCounting_;
228 };
229
230 #endif