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