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