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 * Declares the GridSet class.
41 * This class holds the grids for the local and non-local domain decomposition
42 * zones, as well as the cell and atom data that covers all grids.
44 * \author Berk Hess <hess@kth.se>
45 * \ingroup module_nbnxm
48 #ifndef GMX_NBNXM_GRIDSET_H
49 #define GMX_NBNXM_GRIDSET_H
54 #include "gromacs/math/vec.h"
55 #include "gromacs/math/vectypes.h"
56 #include "gromacs/utility/alignedallocator.h"
57 #include "gromacs/utility/arrayref.h"
60 #include "gridsetdata.h"
63 struct nbnxn_atomdata_t
;
64 enum class PairlistType
;
68 class UpdateGroupsCog
;
75 * \brief Holds a set of search grids for the local + non-local DD zones
81 * \brief Description of the domain setup: PBC and the connections between domains
85 //! Constructor, without DD \p numDDCells and \p ddZones should be nullptr
87 const ivec
*numDDCells
,
88 const gmx_domdec_zones_t
*ddZones
);
92 //! Are there multiple domains?
93 bool haveMultipleDomains
;
94 //! Are there multiple domains along each dimension?
95 std::array
<bool, DIM
> haveMultipleDomainsPerDim
;
96 //! The domain decomposition zone setup
97 const gmx_domdec_zones_t
*zones
;
100 //! Constructs a grid set for 1 or multiple DD zones, when numDDCells!=nullptr
102 const ivec
*numDDCells
,
103 const gmx_domdec_zones_t
*ddZones
,
104 PairlistType pairlistType
,
107 gmx::PinningPolicy pinningPolicy
);
109 //! Puts the atoms in \p ddZone on the grid and copies the coordinates to \p nbat
110 void putOnGrid(const matrix box
,
112 const rvec lowerCorner
,
113 const rvec upperCorner
,
114 const gmx::UpdateGroupsCog
*updateGroupsCog
,
118 gmx::ArrayRef
<const int> atomInfo
,
119 gmx::ArrayRef
<const gmx::RVec
> x
,
122 nbnxn_atomdata_t
*nbat
);
124 //! Returns the domain setup
125 const DomainSetup
domainSetup() const
130 //! Returns the total number of atoms in the grid set, including padding
131 int numGridAtomsTotal() const
133 return grids_
.back().atomIndexEnd();
136 //! Returns the number of local real atoms, i.e. without padded atoms
137 int numRealAtomsLocal() const
139 return numRealAtomsLocal_
;
142 //! Returns the number of total real atoms, i.e. without padded atoms
143 int numRealAtomsTotal() const
145 return numRealAtomsTotal_
;
148 //! Returns the atom order on the grid for the local atoms
149 gmx::ArrayRef
<const int> getLocalAtomorder() const
151 /* Return the atom order for the home cell (index 0) */
152 const int numIndices
= grids_
[0].atomIndexEnd() - grids_
[0].firstAtomInColumn(0);
154 return gmx::constArrayRefFromArray(atomIndices().data(), numIndices
);
157 //! Sets the order of the local atoms to the order grid atom ordering
158 void setLocalAtomOrder();
160 //! Returns the list of grids
161 gmx::ArrayRef
<const Grid
> grids() const
166 //! Returns the grid atom indices covering all grids
167 gmx::ArrayRef
<const int> cells() const
169 return gridSetData_
.cells
;
172 //! Returns the grid atom indices covering all grids
173 gmx::ArrayRef
<const int> atomIndices() const
175 return gridSetData_
.atomIndices
;
178 //! Returns whether we have perturbed non-bonded interactions
184 //! Returns the unit cell in \p box
185 void getBox(matrix box
) const
190 //! Returns the maximum number of columns across all grids
191 int numColumnsMax() const
193 return numColumnsMax_
;
196 //! Sets the maximum number of columns across all grids
197 void setNumColumnsMax(int numColumnsMax
)
199 numColumnsMax_
= numColumnsMax
;
205 DomainSetup domainSetup_
;
207 std::vector
<Grid
> grids_
;
208 //! The cell and atom index data which runs over all grids
209 GridSetData gridSetData_
;
210 //! Tells whether we have perturbed non-bonded interactions
212 //! The periodic unit-cell
214 //! The number of local real atoms, i.e. without padded atoms, local atoms: 0 to numAtomsLocal_
215 int numRealAtomsLocal_
;
216 //! The total number of real atoms, i.e. without padded atoms
217 int numRealAtomsTotal_
;
218 //! Working data for constructing a single grid, one entry per thread
219 std::vector
<GridWork
> gridWork_
;
220 //! Maximum number of columns across all grids