Enable TPI with the Verlet cut-off scheme
[gromacs.git] / src / gromacs / nbnxm / gridset.cpp
blob28a88bbeb90f0560f76f09c850c872d686db8ab9
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 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.
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.
36 /*! \internal \file
38 * \brief
39 * Implements the GridSet class.
41 * \author Berk Hess <hess@kth.se>
42 * \ingroup module_nbnxm
45 #include "gmxpre.h"
47 #include "gridset.h"
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"
54 namespace Nbnxm
57 //! Returns the number of search grids
58 static int numGrids(const GridSet::DomainSetup &domainSetup)
60 int numGrids;
61 if (domainSetup.doTestParticleInsertion)
63 numGrids = 2;
65 else
67 numGrids = 1;
68 for (auto haveDD : domainSetup.haveMultipleDomainsPerDim)
70 if (haveDD)
72 numGrids *= 2;
77 return numGrids;
80 GridSet::DomainSetup::DomainSetup(const int ePBC,
81 const bool doTestParticleInsertion,
82 const ivec *numDDCells,
83 const gmx_domdec_zones_t *ddZones) :
84 ePBC(ePBC),
85 doTestParticleInsertion(doTestParticleInsertion),
86 haveMultipleDomains(numDDCells != nullptr),
87 zones(ddZones)
89 for (int d = 0; d < DIM; d++)
91 haveMultipleDomainsPerDim[d] = (numDDCells != nullptr && (*numDDCells)[d] > 1);
95 GridSet::GridSet(const int ePBC,
96 const bool doTestParticleInsertion,
97 const ivec *numDDCells,
98 const gmx_domdec_zones_t *ddZones,
99 const PairlistType pairlistType,
100 const bool haveFep,
101 const int numThreads,
102 gmx::PinningPolicy pinningPolicy) :
103 domainSetup_(ePBC, doTestParticleInsertion, numDDCells, ddZones),
104 grids_(numGrids(domainSetup_), Grid(pairlistType, haveFep_)),
105 haveFep_(haveFep),
106 numRealAtomsLocal_(0),
107 numRealAtomsTotal_(0),
108 gridWork_(numThreads)
110 clear_mat(box_);
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];
120 int atomIndex = 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;
129 atomIndex++;
130 cellIndex++;
135 // TODO: Move to gridset.cpp
136 void GridSet::putOnGrid(const matrix box,
137 const int gridIndex,
138 const rvec lowerCorner,
139 const rvec upperCorner,
140 const gmx::UpdateGroupsCog *updateGroupsCog,
141 const int atomStart,
142 const int atomEnd,
143 real atomDensity,
144 gmx::ArrayRef<const int> atomInfo,
145 gmx::ArrayRef<const gmx::RVec> x,
146 const int numAtomsMoved,
147 const int *move,
148 nbnxn_atomdata_t *nbat)
150 Nbnxm::Grid &grid = grids_[gridIndex];
152 int cellOffset;
153 if (gridIndex == 0)
155 cellOffset = 0;
157 else
159 const Nbnxm::Grid &previousGrid = grids_[gridIndex - 1];
160 cellOffset = previousGrid.atomIndexEnd()/previousGrid.geometry().numAtomsPerCell;
163 const int n = atomEnd - atomStart;
165 real maxAtomGroupRadius;
166 if (gridIndex == 0)
168 copy_mat(box, box_);
170 numRealAtomsLocal_ = atomEnd - numAtomsMoved;
171 /* We assume that nbnxn_put_on_grid is called first
172 * for the local atoms (gridIndex=0).
174 numRealAtomsTotal_ = atomEnd - numAtomsMoved;
176 maxAtomGroupRadius = (updateGroupsCog ? updateGroupsCog->maxUpdateGroupRadius() : 0);
178 if (debug)
180 fprintf(debug, "natoms_local = %5d atom_density = %5.1f\n",
181 numRealAtomsLocal_, atomDensity);
184 else
186 const Nbnxm::Grid::Dimensions &dimsGrid0 = grids_[0].dimensions();
187 atomDensity = dimsGrid0.atomDensity;
188 maxAtomGroupRadius = dimsGrid0.maxAtomGroupRadius;
190 numRealAtomsTotal_ = std::max(numRealAtomsTotal_, atomEnd);
193 /* We always use the home zone (grid[0]) for setting the cell size,
194 * since determining densities for non-local zones is difficult.
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(ddZone, n - numAtomsMoved,
200 lowerCorner, upperCorner,
201 atomDensity,
202 maxAtomGroupRadius,
203 haveFep_,
204 pinPolicy);
206 for (GridWork &work : gridWork_)
208 work.numAtomsPerColumn.resize(grid.numColumns() + 1);
211 /* Make space for the new cell indices */
212 gridSetData_.cells.resize(atomEnd);
214 const int nthread = gmx_omp_nthreads_get(emntPairsearch);
216 #pragma omp parallel for num_threads(nthread) schedule(static)
217 for (int thread = 0; thread < nthread; thread++)
221 Grid::calcColumnIndices(grid.dimensions(),
222 updateGroupsCog,
223 atomStart, atomEnd, x,
224 ddZone, move, thread, nthread,
225 gridSetData_.cells,
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 */
232 grid.setCellIndices(ddZone, cellOffset, &gridSetData_, gridWork_,
233 atomStart, atomEnd, atomInfo.data(), x, numAtomsMoved, nbat);
235 if (gridIndex == 0)
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 (const auto &grid : grids())
248 maxNumColumns = std::max(maxNumColumns, grid.numColumns());
250 setNumColumnsMax(maxNumColumns);
254 } // namespace Nbnxm