Use gmx::Range in Nbnxm gridding functions
[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, 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/utility/fatalerror.h"
52
53 #include "atomdata.h"
54
55 namespace Nbnxm
56 {
57
58 //! Returns the number of search grids
59 static int numGrids(const GridSet::DomainSetup &domainSetup)
60 {
61     int numGrids;
62     if (domainSetup.doTestParticleInsertion)
63     {
64         numGrids = 2;
65     }
66     else
67     {
68         numGrids = 1;
69         for (auto haveDD : domainSetup.haveMultipleDomainsPerDim)
70         {
71             if (haveDD)
72             {
73                 numGrids *= 2;
74             }
75         }
76     }
77
78     return numGrids;
79 }
80
81 GridSet::DomainSetup::DomainSetup(const int                 ePBC,
82                                   const bool                doTestParticleInsertion,
83                                   const ivec               *numDDCells,
84                                   const gmx_domdec_zones_t *ddZones) :
85     ePBC(ePBC),
86     doTestParticleInsertion(doTestParticleInsertion),
87     haveMultipleDomains(numDDCells != nullptr),
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 int                 ePBC,
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_(ePBC, 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 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 gmx::Range<int>           atomRange,
142                         real                            atomDensity,
143                         gmx::ArrayRef<const int>        atomInfo,
144                         gmx::ArrayRef<const gmx::RVec>  x,
145                         const int                       numAtomsMoved,
146                         const int                      *move,
147                         nbnxn_atomdata_t               *nbat)
148 {
149     Nbnxm::Grid  &grid = grids_[gridIndex];
150
151     int           cellOffset;
152     if (gridIndex == 0)
153     {
154         cellOffset = 0;
155     }
156     else
157     {
158         const Nbnxm::Grid &previousGrid = grids_[gridIndex - 1];
159         cellOffset = previousGrid.atomIndexEnd()/previousGrid.geometry().numAtomsPerCell;
160     }
161
162     const int n = atomRange.size();
163
164     real      maxAtomGroupRadius;
165     if (gridIndex == 0)
166     {
167         copy_mat(box, box_);
168
169         numRealAtomsLocal_ = *atomRange.end() - numAtomsMoved;
170         /* We assume that nbnxn_put_on_grid is called first
171          * for the local atoms (gridIndex=0).
172          */
173         numRealAtomsTotal_ = *atomRange.end() - numAtomsMoved;
174
175         maxAtomGroupRadius = (updateGroupsCog ? updateGroupsCog->maxUpdateGroupRadius() : 0);
176
177         if (debug)
178         {
179             fprintf(debug, "natoms_local = %5d atom_density = %5.1f\n",
180                     numRealAtomsLocal_, atomDensity);
181         }
182     }
183     else
184     {
185         const Nbnxm::Grid::Dimensions &dimsGrid0 = grids_[0].dimensions();
186         atomDensity        = dimsGrid0.atomDensity;
187         maxAtomGroupRadius = dimsGrid0.maxAtomGroupRadius;
188
189         numRealAtomsTotal_ = std::max(numRealAtomsTotal_, *atomRange.end());
190     }
191
192     /* We always use the home zone (grid[0]) for setting the cell size,
193      * since determining densities for non-local zones is difficult.
194      */
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();
198     grid.setDimensions(ddZone, n - numAtomsMoved,
199                        lowerCorner, upperCorner,
200                        atomDensity,
201                        maxAtomGroupRadius,
202                        haveFep_,
203                        pinPolicy);
204
205     for (GridWork &work : gridWork_)
206     {
207         work.numAtomsPerColumn.resize(grid.numColumns() + 1);
208     }
209
210     /* Make space for the new cell indices */
211     gridSetData_.cells.resize(*atomRange.end());
212
213     const int nthread = gmx_omp_nthreads_get(emntPairsearch);
214     GMX_ASSERT(nthread > 0, "We expect the OpenMP thread count to be set");
215
216 #pragma omp parallel for num_threads(nthread) schedule(static)
217     for (int thread = 0; thread < nthread; thread++)
218     {
219         try
220         {
221             Grid::calcColumnIndices(grid.dimensions(),
222                                     updateGroupsCog,
223                                     atomRange, x,
224                                     ddZone, move, thread, nthread,
225                                     gridSetData_.cells,
226                                     gridWork_[thread].numAtomsPerColumn);
227         }
228         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
229     }
230
231     /* Copy the already computed cell indices to the grid and sort, when needed */
232     grid.setCellIndices(ddZone, cellOffset, &gridSetData_, gridWork_,
233                         atomRange, atomInfo.data(), x, numAtomsMoved, nbat);
234
235     if (gridIndex == 0)
236     {
237         nbat->natoms_local = nbat->numAtoms();
238     }
239     if (gridIndex == gmx::ssize(grids_) - 1)
240     {
241         /* We are done setting up all grids, we can resize the force buffers */
242         nbat->resizeForceBuffers();
243     }
244
245     int maxNumColumns = 0;
246     for (int i = 0; i <= gridIndex; i++)
247     {
248         maxNumColumns = std::max(maxNumColumns, grids_[i].numColumns());
249     }
250     setNumColumnsMax(maxNumColumns);
251 }
252
253 } // namespace Nbnxm