Unify handling of GMX_ENABLE_GPU_TIMING and GMX_DISABLE_GPU_TIMING
[alexxy/gromacs.git] / src / gromacs / nbnxm / gridset.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2019,2020,2021, 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
39  * Implements the GridSet class.
40  *
41  * \author Berk Hess <hess@kth.se>
42  * \ingroup module_nbnxm
43  */
44
45 #include "gmxpre.h"
46
47 #include "gridset.h"
48
49 #include "gromacs/mdlib/gmx_omp_nthreads.h"
50 #include "gromacs/mdlib/updategroupscog.h"
51 #include "gromacs/nbnxm/atomdata.h"
52 #include "gromacs/utility/fatalerror.h"
53
54 namespace Nbnxm
55 {
56
57 //! Returns the number of search grids
58 static int numGrids(const GridSet::DomainSetup& domainSetup)
59 {
60     // One grid for the test particle, one for the rest
61     static constexpr int sc_numGridsForTestParticleInsertion = 2;
62     if (domainSetup.doTestParticleInsertion)
63     {
64         return sc_numGridsForTestParticleInsertion;
65     }
66     else
67     {
68         int numGrids = 1;
69         for (auto haveDD : domainSetup.haveMultipleDomainsPerDim)
70         {
71             if (haveDD)
72             {
73                 numGrids *= 2;
74             }
75         }
76         return numGrids;
77     }
78 }
79
80 GridSet::DomainSetup::DomainSetup(const PbcType             pbcType,
81                                   const bool                doTestParticleInsertion,
82                                   const ivec*               numDDCells,
83                                   const gmx_domdec_zones_t* ddZones) :
84     pbcType(pbcType),
85     doTestParticleInsertion(doTestParticleInsertion),
86     haveMultipleDomains(numDDCells != nullptr
87                         && (*numDDCells)[XX] * (*numDDCells)[YY] * (*numDDCells)[ZZ] > 1),
88     zones(ddZones)
89 {
90     for (int d = 0; d < DIM; d++)
91     {
92         haveMultipleDomainsPerDim[d] = (numDDCells != nullptr && (*numDDCells)[d] > 1);
93     }
94 }
95
96 GridSet::GridSet(const PbcType             pbcType,
97                  const bool                doTestParticleInsertion,
98                  const ivec*               numDDCells,
99                  const gmx_domdec_zones_t* ddZones,
100                  const PairlistType        pairlistType,
101                  const bool                haveFep,
102                  const int                 numThreads,
103                  gmx::PinningPolicy        pinningPolicy) :
104     domainSetup_(pbcType, doTestParticleInsertion, numDDCells, ddZones),
105     grids_(numGrids(domainSetup_), Grid(pairlistType, haveFep_)),
106     haveFep_(haveFep),
107     numRealAtomsLocal_(0),
108     numRealAtomsTotal_(0),
109     gridWork_(numThreads)
110 {
111     clear_mat(box_);
112     changePinningPolicy(&gridSetData_.cells, pinningPolicy);
113     changePinningPolicy(&gridSetData_.atomIndices, pinningPolicy);
114 }
115
116 void GridSet::setLocalAtomOrder()
117 {
118     /* Set the atom order for the home cell (index 0) */
119     const Nbnxm::Grid& grid = grids_[0];
120
121     int atomIndex = 0;
122     for (int cxy = 0; cxy < grid.numColumns(); cxy++)
123     {
124         const int numAtoms  = grid.numAtomsInColumn(cxy);
125         int       cellIndex = grid.firstCellInColumn(cxy) * grid.geometry().numAtomsPerCell;
126         for (int i = 0; i < numAtoms; i++)
127         {
128             gridSetData_.atomIndices[cellIndex] = atomIndex;
129             gridSetData_.cells[atomIndex]       = cellIndex;
130             atomIndex++;
131             cellIndex++;
132         }
133     }
134 }
135
136 static int getGridOffset(gmx::ArrayRef<const Grid> grids, int gridIndex)
137 {
138     if (gridIndex == 0)
139     {
140         return 0;
141     }
142     else
143     {
144         const Nbnxm::Grid& previousGrid = grids[gridIndex - 1];
145         return previousGrid.atomIndexEnd() / previousGrid.geometry().numAtomsPerCell;
146     }
147 }
148
149 void GridSet::putOnGrid(const matrix                   box,
150                         const int                      gridIndex,
151                         const rvec                     lowerCorner,
152                         const rvec                     upperCorner,
153                         const gmx::UpdateGroupsCog*    updateGroupsCog,
154                         const gmx::Range<int>          atomRange,
155                         real                           atomDensity,
156                         gmx::ArrayRef<const int64_t>   atomInfo,
157                         gmx::ArrayRef<const gmx::RVec> x,
158                         const int                      numAtomsMoved,
159                         const int*                     move,
160                         nbnxn_atomdata_t*              nbat)
161 {
162     Nbnxm::Grid& grid               = grids_[gridIndex];
163     const int    cellOffset         = getGridOffset(grids_, gridIndex);
164     const int    n                  = atomRange.size();
165     real         maxAtomGroupRadius = NAN;
166
167     if (gridIndex == 0)
168     {
169         copy_mat(box, box_);
170
171         numRealAtomsLocal_ = *atomRange.end() - numAtomsMoved;
172         /* We assume that nbnxn_put_on_grid is called first
173          * for the local atoms (gridIndex=0).
174          */
175         numRealAtomsTotal_ = *atomRange.end() - numAtomsMoved;
176
177         maxAtomGroupRadius = (updateGroupsCog ? updateGroupsCog->maxUpdateGroupRadius() : 0);
178
179         if (debug)
180         {
181             fprintf(debug, "natoms_local = %5d atom_density = %5.1f\n", numRealAtomsLocal_, atomDensity);
182         }
183     }
184     else
185     {
186         const Nbnxm::Grid::Dimensions& dimsGrid0 = grids_[0].dimensions();
187         atomDensity                              = dimsGrid0.atomDensity;
188         maxAtomGroupRadius                       = dimsGrid0.maxAtomGroupRadius;
189
190         numRealAtomsTotal_ = std::max(numRealAtomsTotal_, *atomRange.end());
191     }
192
193     /* We always use the home zone (grid[0]) for setting the cell size,
194      * since determining densities for non-local zones is difficult.
195      */
196     const int ddZone = (domainSetup_.doTestParticleInsertion ? 0 : gridIndex);
197     // grid data used in GPU transfers inherits the gridset pinning policy
198     auto pinPolicy = gridSetData_.cells.get_allocator().pinningPolicy();
199     grid.setDimensions(
200             ddZone, n - numAtomsMoved, lowerCorner, upperCorner, atomDensity, maxAtomGroupRadius, haveFep_, pinPolicy);
201
202     for (GridWork& work : gridWork_)
203     {
204         work.numAtomsPerColumn.resize(grid.numColumns() + 1);
205     }
206
207     /* Make space for the new cell indices */
208     gridSetData_.cells.resize(*atomRange.end());
209
210     const int nthread = gmx_omp_nthreads_get(ModuleMultiThread::Pairsearch);
211     GMX_ASSERT(nthread > 0, "We expect the OpenMP thread count to be set");
212
213 #pragma omp parallel for num_threads(nthread) schedule(static)
214     for (int thread = 0; thread < nthread; thread++)
215     {
216         try
217         {
218             Grid::calcColumnIndices(grid.dimensions(),
219                                     updateGroupsCog,
220                                     atomRange,
221                                     x,
222                                     ddZone,
223                                     move,
224                                     thread,
225                                     nthread,
226                                     gridSetData_.cells,
227                                     gridWork_[thread].numAtomsPerColumn);
228         }
229         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
230     }
231
232     /* Copy the already computed cell indices to the grid and sort, when needed */
233     grid.setCellIndices(
234             ddZone, cellOffset, &gridSetData_, gridWork_, atomRange, atomInfo, x, numAtomsMoved, nbat);
235
236     if (gridIndex == 0)
237     {
238         nbat->natoms_local = nbat->numAtoms();
239     }
240     if (gridIndex == gmx::ssize(grids_) - 1)
241     {
242         /* We are done setting up all grids, we can resize the force buffers */
243         nbat->resizeForceBuffers();
244     }
245
246     int maxNumColumns = 0;
247     for (int i = 0; i <= gridIndex; i++)
248     {
249         maxNumColumns = std::max(maxNumColumns, grids_[i].numColumns());
250     }
251     setNumColumnsMax(maxNumColumns);
252 }
253
254 } // namespace Nbnxm