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.
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 DD zones
58 static int numDDZones(const std::array
<bool, DIM
> &haveMultipleDomainsPerDim
)
61 for (auto haveDD
: haveMultipleDomainsPerDim
)
72 GridSet::DomainSetup::DomainSetup(const int ePBC
,
73 const ivec
*numDDCells
,
74 const gmx_domdec_zones_t
*ddZones
) :
76 haveMultipleDomains(numDDCells
!= nullptr),
79 for (int d
= 0; d
< DIM
; d
++)
81 haveMultipleDomainsPerDim
[d
] = (numDDCells
!= nullptr && (*numDDCells
)[d
] > 1);
85 GridSet::GridSet(const int ePBC
,
86 const ivec
*numDDCells
,
87 const gmx_domdec_zones_t
*ddZones
,
88 const PairlistType pairlistType
,
91 gmx::PinningPolicy pinningPolicy
) :
92 domainSetup_(ePBC
, numDDCells
, ddZones
),
93 grids_(numDDZones(domainSetup_
.haveMultipleDomainsPerDim
), Grid(pairlistType
, haveFep_
)),
95 numRealAtomsLocal_(0),
96 numRealAtomsTotal_(0),
100 changePinningPolicy(&gridSetData_
.cells
, pinningPolicy
);
101 changePinningPolicy(&gridSetData_
.atomIndices
, pinningPolicy
);
104 void GridSet::setLocalAtomOrder()
106 /* Set the atom order for the home cell (index 0) */
107 const Nbnxm::Grid
&grid
= grids_
[0];
110 for (int cxy
= 0; cxy
< grid
.numColumns(); cxy
++)
112 const int numAtoms
= grid
.numAtomsInColumn(cxy
);
113 int cellIndex
= grid
.firstCellInColumn(cxy
)*grid
.geometry().numAtomsPerCell
;
114 for (int i
= 0; i
< numAtoms
; i
++)
116 gridSetData_
.atomIndices
[cellIndex
] = atomIndex
;
117 gridSetData_
.cells
[atomIndex
] = cellIndex
;
124 // TODO: Move to gridset.cpp
125 void GridSet::putOnGrid(const matrix box
,
127 const rvec lowerCorner
,
128 const rvec upperCorner
,
129 const gmx::UpdateGroupsCog
*updateGroupsCog
,
133 gmx::ArrayRef
<const int> atomInfo
,
134 gmx::ArrayRef
<const gmx::RVec
> x
,
135 const int numAtomsMoved
,
137 nbnxn_atomdata_t
*nbat
)
139 Nbnxm::Grid
&grid
= grids_
[ddZone
];
148 const Nbnxm::Grid
&previousGrid
= grids_
[ddZone
- 1];
149 cellOffset
= previousGrid
.atomIndexEnd()/previousGrid
.geometry().numAtomsPerCell
;
152 const int n
= atomEnd
- atomStart
;
154 real maxAtomGroupRadius
;
159 numRealAtomsLocal_
= atomEnd
- numAtomsMoved
;
160 /* We assume that nbnxn_put_on_grid is called first
161 * for the local atoms (ddZone=0).
163 numRealAtomsTotal_
= atomEnd
- numAtomsMoved
;
165 maxAtomGroupRadius
= (updateGroupsCog
? updateGroupsCog
->maxUpdateGroupRadius() : 0);
169 fprintf(debug
, "natoms_local = %5d atom_density = %5.1f\n",
170 numRealAtomsLocal_
, atomDensity
);
175 const Nbnxm::Grid::Dimensions
&dimsGrid0
= grids_
[0].dimensions();
176 atomDensity
= dimsGrid0
.atomDensity
;
177 maxAtomGroupRadius
= dimsGrid0
.maxAtomGroupRadius
;
179 numRealAtomsTotal_
= std::max(numRealAtomsTotal_
, atomEnd
);
182 /* We always use the home zone (grid[0]) for setting the cell size,
183 * since determining densities for non-local zones is difficult.
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
,
194 for (GridWork
&work
: gridWork_
)
196 work
.numAtomsPerColumn
.resize(grid
.numColumns() + 1);
199 /* Make space for the new cell indices */
200 gridSetData_
.cells
.resize(atomEnd
);
202 const int nthread
= gmx_omp_nthreads_get(emntPairsearch
);
204 #pragma omp parallel for num_threads(nthread) schedule(static)
205 for (int thread
= 0; thread
< nthread
; thread
++)
209 Grid::calcColumnIndices(grid
.dimensions(),
211 atomStart
, atomEnd
, x
,
212 ddZone
, move
, thread
, nthread
,
214 gridWork_
[thread
].numAtomsPerColumn
);
216 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
;
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
);
225 nbat
->natoms_local
= nbat
->numAtoms();
227 if (ddZone
== gmx::ssize(grids_
) - 1)
229 /* We are done setting up all grids, we can resize the force buffers */
230 nbat
->resizeForceBuffers();
233 int maxNumColumns
= 0;
234 for (const auto &grid
: grids())
236 maxNumColumns
= std::max(maxNumColumns
, grid
.numColumns());
238 setNumColumnsMax(maxNumColumns
);