34369acf9498be2d1dbc67a48765827d3c1a2211
[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/nbnxm/atomdata.h"
52 #include "gromacs/utility/fatalerror.h"
53
54 namespace Nbnxm
55 {
56
57 //! Returns the number of DD zones
58 static int numDDZones(const std::array<bool, DIM> &haveMultipleDomainsPerDim)
59 {
60     int  numDDZones = 1;
61     for (auto haveDD : haveMultipleDomainsPerDim)
62     {
63         if (haveDD)
64         {
65             numDDZones *= 2;
66         }
67     }
68
69     return numDDZones;
70 }
71
72 GridSet::DomainSetup::DomainSetup(const int                 ePBC,
73                                   const ivec               *numDDCells,
74                                   const gmx_domdec_zones_t *ddZones) :
75     ePBC(ePBC),
76     haveMultipleDomains(numDDCells != nullptr),
77     zones(ddZones)
78 {
79     for (int d = 0; d < DIM; d++)
80     {
81         haveMultipleDomainsPerDim[d] = (numDDCells != nullptr && (*numDDCells)[d] > 1);
82     }
83 }
84
85 GridSet::GridSet(const int                 ePBC,
86                  const ivec               *numDDCells,
87                  const gmx_domdec_zones_t *ddZones,
88                  const PairlistType        pairlistType,
89                  const bool                haveFep,
90                  const int                 numThreads,
91                  gmx::PinningPolicy        pinningPolicy) :
92     domainSetup_(ePBC, numDDCells, ddZones),
93     grids_(numDDZones(domainSetup_.haveMultipleDomainsPerDim), Grid(pairlistType, haveFep_)),
94     haveFep_(haveFep),
95     numRealAtomsLocal_(0),
96     numRealAtomsTotal_(0),
97     gridWork_(numThreads)
98 {
99     clear_mat(box_);
100     changePinningPolicy(&gridSetData_.cells, pinningPolicy);
101     changePinningPolicy(&gridSetData_.atomIndices, pinningPolicy);
102 }
103
104 void GridSet::setLocalAtomOrder()
105 {
106     /* Set the atom order for the home cell (index 0) */
107     const Nbnxm::Grid &grid = grids_[0];
108
109     int                atomIndex = 0;
110     for (int cxy = 0; cxy < grid.numColumns(); cxy++)
111     {
112         const int numAtoms  = grid.numAtomsInColumn(cxy);
113         int       cellIndex = grid.firstCellInColumn(cxy)*grid.geometry().numAtomsPerCell;
114         for (int i = 0; i < numAtoms; i++)
115         {
116             gridSetData_.atomIndices[cellIndex] = atomIndex;
117             gridSetData_.cells[atomIndex]       = cellIndex;
118             atomIndex++;
119             cellIndex++;
120         }
121     }
122 }
123
124 // TODO: Move to gridset.cpp
125 void GridSet::putOnGrid(const matrix                    box,
126                         const int                       ddZone,
127                         const rvec                      lowerCorner,
128                         const rvec                      upperCorner,
129                         const gmx::UpdateGroupsCog     *updateGroupsCog,
130                         const int                       atomStart,
131                         const int                       atomEnd,
132                         real                            atomDensity,
133                         gmx::ArrayRef<const int>        atomInfo,
134                         gmx::ArrayRef<const gmx::RVec>  x,
135                         const int                       numAtomsMoved,
136                         const int                      *move,
137                         nbnxn_atomdata_t               *nbat)
138 {
139     Nbnxm::Grid  &grid = grids_[ddZone];
140
141     int           cellOffset;
142     if (ddZone == 0)
143     {
144         cellOffset = 0;
145     }
146     else
147     {
148         const Nbnxm::Grid &previousGrid = grids_[ddZone - 1];
149         cellOffset = previousGrid.atomIndexEnd()/previousGrid.geometry().numAtomsPerCell;
150     }
151
152     const int n = atomEnd - atomStart;
153
154     real      maxAtomGroupRadius;
155     if (ddZone == 0)
156     {
157         copy_mat(box, box_);
158
159         numRealAtomsLocal_ = atomEnd - numAtomsMoved;
160         /* We assume that nbnxn_put_on_grid is called first
161          * for the local atoms (ddZone=0).
162          */
163         numRealAtomsTotal_ = atomEnd - numAtomsMoved;
164
165         maxAtomGroupRadius = (updateGroupsCog ? updateGroupsCog->maxUpdateGroupRadius() : 0);
166
167         if (debug)
168         {
169             fprintf(debug, "natoms_local = %5d atom_density = %5.1f\n",
170                     numRealAtomsLocal_, atomDensity);
171         }
172     }
173     else
174     {
175         const Nbnxm::Grid::Dimensions &dimsGrid0 = grids_[0].dimensions();
176         atomDensity        = dimsGrid0.atomDensity;
177         maxAtomGroupRadius = dimsGrid0.maxAtomGroupRadius;
178
179         numRealAtomsTotal_ = std::max(numRealAtomsTotal_, atomEnd);
180     }
181
182     /* We always use the home zone (grid[0]) for setting the cell size,
183      * since determining densities for non-local zones is difficult.
184      */
185     // grid data used in GPU transfers inherits the gridset pinnin policy
186     auto pinPolicy = gridSetData_.cells.get_allocator().pinningPolicy();
187     grid.setDimensions(ddZone, n - numAtomsMoved,
188                        lowerCorner, upperCorner,
189                        atomDensity,
190                        maxAtomGroupRadius,
191                        haveFep_,
192                        pinPolicy);
193
194     for (GridWork &work : gridWork_)
195     {
196         work.numAtomsPerColumn.resize(grid.numColumns() + 1);
197     }
198
199     /* Make space for the new cell indices */
200     gridSetData_.cells.resize(atomEnd);
201
202     const int nthread = gmx_omp_nthreads_get(emntPairsearch);
203
204 #pragma omp parallel for num_threads(nthread) schedule(static)
205     for (int thread = 0; thread < nthread; thread++)
206     {
207         try
208         {
209             Grid::calcColumnIndices(grid.dimensions(),
210                                     updateGroupsCog,
211                                     atomStart, atomEnd, x,
212                                     ddZone, move, thread, nthread,
213                                     gridSetData_.cells,
214                                     gridWork_[thread].numAtomsPerColumn);
215         }
216         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
217     }
218
219     /* Copy the already computed cell indices to the grid and sort, when needed */
220     grid.setCellIndices(ddZone, cellOffset, &gridSetData_, gridWork_,
221                         atomStart, atomEnd, atomInfo.data(), x, numAtomsMoved, nbat);
222
223     if (ddZone == 0)
224     {
225         nbat->natoms_local = nbat->numAtoms();
226     }
227     if (ddZone == gmx::ssize(grids_) - 1)
228     {
229         /* We are done setting up all grids, we can resize the force buffers */
230         nbat->resizeForceBuffers();
231     }
232
233     int maxNumColumns = 0;
234     for (const auto &grid : grids())
235     {
236         maxNumColumns = std::max(maxNumColumns, grid.numColumns());
237     }
238     setNumColumnsMax(maxNumColumns);
239
240 }
241
242 } // namespace Nbnxm