Clean up mdrun help description
[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 // TODO: Move nbnxn_search_work_t definition to its own file
114
115 /* Thread-local work struct, contains working data for Grid */
116 struct PairsearchWork
117 {
118     PairsearchWork();
119
120     ~PairsearchWork();
121
122     gmx_cache_protect_t       cp0;          /* Buffer to avoid cache polution */
123
124     std::vector<int>          sortBuffer;   /* Temporary buffer for sorting atoms within a grid column */
125
126     nbnxn_buffer_flags_t      buffer_flags; /* Flags for force buffer access */
127
128     int                       ndistc;       /* Number of distance checks for flop counting */
129
130
131     std::unique_ptr<t_nblist> nbl_fep;      /* Temporary FEP list for load balancing */
132
133     nbnxn_cycle_t             cycleCounter; /* Counter for thread-local cycles */
134
135     gmx_cache_protect_t       cp1;          /* Buffer to avoid cache polution */
136 };
137
138 /* Main pair-search struct, contains the grid(s), not the pair-list(s) */
139 class PairSearch
140 {
141     public:
142         /*! \internal
143          * \brief Description of the domain setup: PBC and the connections between domains
144          */
145         struct DomainSetup
146         {
147             /*! \internal
148              * \brief Description of the domain setup: PBC and the connections between domains
149              */
150             //! Constructor, without DD \p numDDCells and \p ddZones should be nullptr
151             DomainSetup(int                       ePBC,
152                         const ivec               *numDDCells,
153                         const gmx_domdec_zones_t *ddZones);
154
155             //! The type of PBC
156             int                       ePBC;
157             //! Tells whether we are using domain decomposition
158             bool                      haveDomDec;
159             //! Tells whether we are using domain decomposition per dimension
160             std::array<bool, DIM>     haveDomDecPerDim;
161             //! The domain decomposition zone setup
162             const gmx_domdec_zones_t *zones;
163         };
164
165         //! Local cycle count enum for profiling different parts of search
166         enum {
167             enbsCCgrid, enbsCCsearch, enbsCCcombine, enbsCCnr
168         };
169
170         struct SearchCycleCounting
171         {
172             //! Start a pair search cycle counter
173             void start(const int enbsCC)
174             {
175                 cc_[enbsCC].start();
176             }
177
178             //! Stop a pair search cycle counter
179             void stop(const int enbsCC)
180             {
181                 cc_[enbsCC].stop();
182             }
183
184             //! Print the cycle counts to \p fp
185             void printCycles(FILE                               *fp,
186                              gmx::ArrayRef<const PairsearchWork> work) const;
187
188             bool          recordCycles_ = false;
189             int           searchCount_  = 0;
190             nbnxn_cycle_t cc_[enbsCCnr];
191         };
192
193         //! Puts the atoms in \p ddZone on the grid and copies the coordinates to \p nbat
194         void putOnGrid(const matrix                    box,
195                        int                             ddZone,
196                        const rvec                      lowerCorner,
197                        const rvec                      upperCorner,
198                        const gmx::UpdateGroupsCog     *updateGroupsCog,
199                        int                             atomStart,
200                        int                             atomEnd,
201                        real                            atomDensity,
202                        const int                      *atinfo,
203                        gmx::ArrayRef<const gmx::RVec>  x,
204                        int                             numAtomsMoved,
205                        const int                      *move,
206                        nbnxn_atomdata_t               *nbat)
207         {
208             cycleCounting_.start(enbsCCgrid);
209
210             gridSet_.putOnGrid(box, ddZone, lowerCorner, upperCorner,
211                                updateGroupsCog, atomStart, atomEnd, atomDensity,
212                                atinfo, x, numAtomsMoved, move, nbat);
213
214             cycleCounting_.stop(enbsCCgrid);
215         }
216
217         /* \brief Constructor
218          *
219          * \param[in] ePBC            The periodic boundary conditions
220          * \param[in] numDDCells      The number of domain decomposition cells per dimension, without DD nullptr should be passed
221          * \param[in] zones           The domain decomposition zone setup, without DD nullptr should be passed
222          * \param[in] haveFep         Tells whether non-bonded interactions are perturbed
223          * \param[in] maxNumThreads   The maximum number of threads used in the search
224          */
225         PairSearch(int                       ePBC,
226                    const ivec               *numDDCells,
227                    const gmx_domdec_zones_t *zones,
228                    PairlistType              pairlistType,
229                    bool                      haveFep,
230                    int                       maxNumthreads);
231
232         //! Sets the order of the local atoms to the order grid atom ordering
233         void setLocalAtomOrder()
234         {
235             gridSet_.setLocalAtomOrder();
236         }
237
238         const DomainSetup domainSetup() const
239         {
240             return domainSetup_;
241         }
242
243         //! Returns the set of search grids
244         const Nbnxm::GridSet &gridSet() const
245         {
246             return gridSet_;
247         }
248
249         //! Returns the list of thread-local work objects
250         gmx::ArrayRef<const PairsearchWork> work() const
251         {
252             return work_;
253         }
254
255         //! Returns the list of thread-local work objects
256         gmx::ArrayRef<PairsearchWork> work()
257         {
258             return work_;
259         }
260
261     private:
262         //! The domain setup
263         DomainSetup                 domainSetup_;
264         //! The set of search grids
265         Nbnxm::GridSet              gridSet_;
266         //! Work objects, one entry for each thread
267         std::vector<PairsearchWork> work_;
268
269     public:
270         //! Cycle counting for measuring components of the search
271         SearchCycleCounting  cycleCounting_;
272 };
273
274 #endif