Make nbnxm a proper module
[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()
81         {
82             start_ = gmx_cycles_read();
83         }
84
85         void stop()
86         {
87             cycles_ += gmx_cycles_read() - start_;
88             count_++;
89         }
90
91         int count() const
92         {
93             return count_;
94         }
95
96         double averageMCycles() const
97         {
98             if (count_ > 0)
99             {
100                 return static_cast<double>(cycles_)*1e-6/count_;
101             }
102             else
103             {
104                 return 0;
105             }
106         }
107
108     private:
109         int          count_  = 0;
110         gmx_cycles_t cycles_ = 0;
111         gmx_cycles_t start_  = 0;
112 };
113
114 //! Local cycle count enum for profiling different parts of search
115 enum {
116     enbsCCgrid, enbsCCsearch, enbsCCcombine, enbsCCnr
117 };
118
119 /*! \internal
120  * \brief Struct for collecting detailed cycle counts for the search
121  */
122 struct SearchCycleCounting
123 {
124     //! Start a pair search cycle counter
125     void start(const int enbsCC)
126     {
127         cc_[enbsCC].start();
128     }
129
130     //! Stop a pair search cycle counter
131     void stop(const int enbsCC)
132     {
133         cc_[enbsCC].stop();
134     }
135
136     //! Print the cycle counts to \p fp
137     void printCycles(FILE                               *fp,
138                      gmx::ArrayRef<const PairsearchWork> work) const;
139
140     //! Tells whether we record cycles
141     bool          recordCycles_ = false;
142     //! The number of times pairsearching has been performed, local+non-local count as 1
143     int           searchCount_  = 0;
144     //! The set of cycle counters
145     nbnxn_cycle_t cc_[enbsCCnr];
146 };
147
148 // TODO: Move nbnxn_search_work_t definition to its own file
149
150 /* Thread-local work struct, contains working data for Grid */
151 struct PairsearchWork
152 {
153     PairsearchWork();
154
155     ~PairsearchWork();
156
157     gmx_cache_protect_t       cp0;          /* Buffer to avoid cache polution */
158
159     std::vector<int>          sortBuffer;   /* Temporary buffer for sorting atoms within a grid column */
160
161     nbnxn_buffer_flags_t      buffer_flags; /* Flags for force buffer access */
162
163     int                       ndistc;       /* Number of distance checks for flop counting */
164
165
166     std::unique_ptr<t_nblist> nbl_fep;      /* Temporary FEP list for load balancing */
167
168     nbnxn_cycle_t             cycleCounter; /* Counter for thread-local cycles */
169
170     gmx_cache_protect_t       cp1;          /* Buffer to avoid cache polution */
171 };
172
173 /* Main pair-search struct, contains the grid(s), not the pair-list(s) */
174 class PairSearch
175 {
176     public:
177         //! Puts the atoms in \p ddZone on the grid and copies the coordinates to \p nbat
178         void putOnGrid(const matrix                    box,
179                        int                             ddZone,
180                        const rvec                      lowerCorner,
181                        const rvec                      upperCorner,
182                        const gmx::UpdateGroupsCog     *updateGroupsCog,
183                        int                             atomStart,
184                        int                             atomEnd,
185                        real                            atomDensity,
186                        gmx::ArrayRef<const int>        atomInfo,
187                        gmx::ArrayRef<const gmx::RVec>  x,
188                        int                             numAtomsMoved,
189                        const int                      *move,
190                        nbnxn_atomdata_t               *nbat)
191         {
192             cycleCounting_.start(enbsCCgrid);
193
194             gridSet_.putOnGrid(box, ddZone, lowerCorner, upperCorner,
195                                updateGroupsCog, atomStart, atomEnd, atomDensity,
196                                atomInfo, x, numAtomsMoved, move, nbat);
197
198             cycleCounting_.stop(enbsCCgrid);
199         }
200
201         /* \brief Constructor
202          *
203          * \param[in] ePBC            The periodic boundary conditions
204          * \param[in] numDDCells      The number of domain decomposition cells per dimension, without DD nullptr should be passed
205          * \param[in] zones           The domain decomposition zone setup, without DD nullptr should be passed
206          * \param[in] haveFep         Tells whether non-bonded interactions are perturbed
207          * \param[in] maxNumThreads   The maximum number of threads used in the search
208          */
209         PairSearch(int                       ePBC,
210                    bool                      doTestParticleInsertion,
211                    const ivec               *numDDCells,
212                    const gmx_domdec_zones_t *zones,
213                    PairlistType              pairlistType,
214                    bool                      haveFep,
215                    int                       maxNumthreads,
216                    gmx::PinningPolicy        pinningPolicy);
217
218         //! Sets the order of the local atoms to the order grid atom ordering
219         void setLocalAtomOrder()
220         {
221             gridSet_.setLocalAtomOrder();
222         }
223
224         //! Returns the set of search grids
225         const Nbnxm::GridSet &gridSet() const
226         {
227             return gridSet_;
228         }
229
230         //! Returns the list of thread-local work objects
231         gmx::ArrayRef<const PairsearchWork> work() const
232         {
233             return work_;
234         }
235
236         //! Returns the list of thread-local work objects
237         gmx::ArrayRef<PairsearchWork> work()
238         {
239             return work_;
240         }
241
242     private:
243         //! The set of search grids
244         Nbnxm::GridSet              gridSet_;
245         //! Work objects, one entry for each thread
246         std::vector<PairsearchWork> work_;
247
248     public:
249         //! Cycle counting for measuring components of the search
250         SearchCycleCounting         cycleCounting_;
251 };
252
253 #endif