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"
62 struct nbnxn_atomdata_t;
63 enum class PairlistType;
67 class UpdateGroupsCog;
74 * \brief Holds a set of search grids for the local + non-local DD zones
80 * \brief Description of the domain setup: PBC and the connections between domains
84 //! Constructor, without DD \p numDDCells and \p ddZones should be nullptr
86 const ivec *numDDCells,
87 const gmx_domdec_zones_t *ddZones);
91 //! Are there multiple domains?
92 bool haveMultipleDomains;
93 //! Are there multiple domains along each dimension?
94 std::array<bool, DIM> haveMultipleDomainsPerDim;
95 //! The domain decomposition zone setup
96 const gmx_domdec_zones_t *zones;
99 //! Constructs a grid set for 1 or multiple DD zones, when numDDCells!=nullptr
101 const ivec *numDDCells,
102 const gmx_domdec_zones_t *ddZones,
103 PairlistType pairlistType,
107 //! Puts the atoms in \p ddZone on the grid and copies the coordinates to \p nbat
108 void putOnGrid(const matrix box,
110 const rvec lowerCorner,
111 const rvec upperCorner,
112 const gmx::UpdateGroupsCog *updateGroupsCog,
117 gmx::ArrayRef<const gmx::RVec> x,
120 nbnxn_atomdata_t *nbat);
122 //! Returns the domain setup
123 const DomainSetup domainSetup() const
128 //! Returns the number of cells along x and y for the local grid
129 void getLocalNumCells(int *numCellsX,
130 int *numCellsY) const
132 *numCellsX = grids_[0].dimensions().numCells[XX];
133 *numCellsY = grids_[0].dimensions().numCells[YY];
136 //! Returns the total number of atoms in the grid set, including padding
137 int numGridAtomsTotal() const
139 return grids_.back().atomIndexEnd();
142 //! Returns the number of local real atoms, i.e. without padded atoms
143 int numRealAtomsLocal() const
145 return numRealAtomsLocal_;
148 //! Returns the number of total real atoms, i.e. without padded atoms
149 int numRealAtomsTotal() const
151 return numRealAtomsTotal_;
154 //! Returns the atom order on the grid for the local atoms
155 gmx::ArrayRef<const int> getLocalAtomorder() const
157 /* Return the atom order for the home cell (index 0) */
158 const int numIndices = grids_[0].atomIndexEnd() - grids_[0].firstAtomInColumn(0);
160 return gmx::constArrayRefFromArray(atomIndices_.data(), numIndices);
163 //! Sets the order of the local atoms to the order grid atom ordering
164 void setLocalAtomOrder();
166 //! Returns the list of grids
167 gmx::ArrayRef<const Grid> grids() const
172 //! Returns the grid atom indices covering all grids
173 gmx::ArrayRef<const int> cells() const
178 //! Returns the grid atom indices covering all grids
179 gmx::ArrayRef<const int> atomIndices() const
184 //! Returns whether we have perturbed non-bonded interactions
190 //! Returns the unit cell in \p box
191 void getBox(matrix box) const
197 //! Returns collection of the data that covers all grids
198 const GridSetData getGridSetData()
200 GridSetData gridSetData = { cells_, atomIndices_, haveFep_ };
207 DomainSetup domainSetup_;
209 std::vector<Grid> grids_;
210 //! The actual cell indices for all atoms, covering all grids
211 std::vector<int> cells_;
212 //! The actual array of atom indices, covering all grids
213 std::vector<int> atomIndices_;
214 //! Tells whether we have perturbed non-bonded interactions
216 //! The periodic unit-cell
218 //! The number of local real atoms, i.e. without padded atoms, local atoms: 0 to numAtomsLocal_
219 int numRealAtomsLocal_;
220 //! The total number of real atoms, i.e. without padded atoms
221 int numRealAtomsTotal_;
222 //! Working data for constructing a single grid, one entry per thread
223 std::vector<GridWork> gridWork_;