2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
39 * Implements the GridSet class.
41 * \author Berk Hess <hess@kth.se>
42 * \ingroup module_nbnxm
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"
57 //! Returns the number of search grids
58 static int numGrids(const GridSet::DomainSetup& domainSetup)
60 // One grid for the test particle, one for the rest
61 static constexpr int sc_numGridsForTestParticleInsertion = 2;
62 if (domainSetup.doTestParticleInsertion)
64 return sc_numGridsForTestParticleInsertion;
69 for (auto haveDD : domainSetup.haveMultipleDomainsPerDim)
80 GridSet::DomainSetup::DomainSetup(const PbcType pbcType,
81 const bool doTestParticleInsertion,
82 const ivec* numDDCells,
83 const gmx_domdec_zones_t* ddZones) :
85 doTestParticleInsertion(doTestParticleInsertion),
86 haveMultipleDomains(numDDCells != nullptr),
89 for (int d = 0; d < DIM; d++)
91 haveMultipleDomainsPerDim[d] = (numDDCells != nullptr && (*numDDCells)[d] > 1);
95 GridSet::GridSet(const PbcType pbcType,
96 const bool doTestParticleInsertion,
97 const ivec* numDDCells,
98 const gmx_domdec_zones_t* ddZones,
99 const PairlistType pairlistType,
101 const int numThreads,
102 gmx::PinningPolicy pinningPolicy) :
103 domainSetup_(pbcType, doTestParticleInsertion, numDDCells, ddZones),
104 grids_(numGrids(domainSetup_), Grid(pairlistType, haveFep_)),
106 numRealAtomsLocal_(0),
107 numRealAtomsTotal_(0),
108 gridWork_(numThreads)
111 changePinningPolicy(&gridSetData_.cells, pinningPolicy);
112 changePinningPolicy(&gridSetData_.atomIndices, pinningPolicy);
115 void GridSet::setLocalAtomOrder()
117 /* Set the atom order for the home cell (index 0) */
118 const Nbnxm::Grid& grid = grids_[0];
121 for (int cxy = 0; cxy < grid.numColumns(); cxy++)
123 const int numAtoms = grid.numAtomsInColumn(cxy);
124 int cellIndex = grid.firstCellInColumn(cxy) * grid.geometry().numAtomsPerCell;
125 for (int i = 0; i < numAtoms; i++)
127 gridSetData_.atomIndices[cellIndex] = atomIndex;
128 gridSetData_.cells[atomIndex] = cellIndex;
135 static int getGridOffset(gmx::ArrayRef<const Grid> grids, int gridIndex)
143 const Nbnxm::Grid& previousGrid = grids[gridIndex - 1];
144 return previousGrid.atomIndexEnd() / previousGrid.geometry().numAtomsPerCell;
148 void GridSet::putOnGrid(const matrix box,
150 const rvec lowerCorner,
151 const rvec upperCorner,
152 const gmx::UpdateGroupsCog* updateGroupsCog,
153 const gmx::Range<int> atomRange,
155 gmx::ArrayRef<const int> atomInfo,
156 gmx::ArrayRef<const gmx::RVec> x,
157 const int numAtomsMoved,
159 nbnxn_atomdata_t* nbat)
161 Nbnxm::Grid& grid = grids_[gridIndex];
162 const int cellOffset = getGridOffset(grids_, gridIndex);
163 const int n = atomRange.size();
164 real maxAtomGroupRadius = NAN;
170 numRealAtomsLocal_ = *atomRange.end() - numAtomsMoved;
171 /* We assume that nbnxn_put_on_grid is called first
172 * for the local atoms (gridIndex=0).
174 numRealAtomsTotal_ = *atomRange.end() - numAtomsMoved;
176 maxAtomGroupRadius = (updateGroupsCog ? updateGroupsCog->maxUpdateGroupRadius() : 0);
180 fprintf(debug, "natoms_local = %5d atom_density = %5.1f\n", numRealAtomsLocal_, atomDensity);
185 const Nbnxm::Grid::Dimensions& dimsGrid0 = grids_[0].dimensions();
186 atomDensity = dimsGrid0.atomDensity;
187 maxAtomGroupRadius = dimsGrid0.maxAtomGroupRadius;
189 numRealAtomsTotal_ = std::max(numRealAtomsTotal_, *atomRange.end());
192 /* We always use the home zone (grid[0]) for setting the cell size,
193 * since determining densities for non-local zones is difficult.
195 const int ddZone = (domainSetup_.doTestParticleInsertion ? 0 : gridIndex);
196 // grid data used in GPU transfers inherits the gridset pinning policy
197 auto pinPolicy = gridSetData_.cells.get_allocator().pinningPolicy();
199 ddZone, n - numAtomsMoved, lowerCorner, upperCorner, atomDensity, maxAtomGroupRadius, haveFep_, pinPolicy);
201 for (GridWork& work : gridWork_)
203 work.numAtomsPerColumn.resize(grid.numColumns() + 1);
206 /* Make space for the new cell indices */
207 gridSetData_.cells.resize(*atomRange.end());
209 const int nthread = gmx_omp_nthreads_get(ModuleMultiThread::Pairsearch);
210 GMX_ASSERT(nthread > 0, "We expect the OpenMP thread count to be set");
212 #pragma omp parallel for num_threads(nthread) schedule(static)
213 for (int thread = 0; thread < nthread; thread++)
217 Grid::calcColumnIndices(grid.dimensions(),
226 gridWork_[thread].numAtomsPerColumn);
228 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
231 /* Copy the already computed cell indices to the grid and sort, when needed */
233 ddZone, cellOffset, &gridSetData_, gridWork_, atomRange, atomInfo.data(), x, numAtomsMoved, nbat);
237 nbat->natoms_local = nbat->numAtoms();
239 if (gridIndex == gmx::ssize(grids_) - 1)
241 /* We are done setting up all grids, we can resize the force buffers */
242 nbat->resizeForceBuffers();
245 int maxNumColumns = 0;
246 for (int i = 0; i <= gridIndex; i++)
248 maxNumColumns = std::max(maxNumColumns, grids_[i].numColumns());
250 setNumColumnsMax(maxNumColumns);