2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2019,2020,2021, 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"
58 #include "gromacs/utility/range.h"
61 #include "gridsetdata.h"
63 struct gmx_domdec_zones_t;
64 struct nbnxn_atomdata_t;
65 enum class PairlistType;
66 enum class PbcType : int;
70 class UpdateGroupsCog;
77 * \brief Holds a set of search grids for the local + non-local DD zones
79 * The are three different possible setups:
80 * - a single grid, this is the standard case without domain decomposition
81 * - one grid for each domain decomposition zone
82 * - with test particle insertion there are two grids, one for the system
83 * to insert in and one for the molecule that is inserted
89 * \brief Description of the domain setup: PBC and the connections between domains
93 //! Constructor, without DD \p numDDCells and \p ddZones should be nullptr
94 DomainSetup(PbcType pbcType,
95 bool doTestParticleInsertion,
96 const ivec* numDDCells,
97 const gmx_domdec_zones_t* ddZones);
101 //! Tells whether we are doing test-particle insertion
102 bool doTestParticleInsertion;
103 //! Are there multiple domains?
104 bool haveMultipleDomains;
105 //! Are there multiple domains along each dimension?
106 std::array<bool, DIM> haveMultipleDomainsPerDim;
107 //! The domain decomposition zone setup
108 const gmx_domdec_zones_t* zones;
111 //! Constructs a grid set for 1 or multiple DD zones, when numDDCells!=nullptr
112 GridSet(PbcType pbcType,
113 bool doTestParticleInsertion,
114 const ivec* numDDCells,
115 const gmx_domdec_zones_t* ddZones,
116 PairlistType pairlistType,
119 gmx::PinningPolicy pinningPolicy);
121 //! Puts the atoms on the grid with index \p gridIndex and copies the coordinates to \p nbat
122 void putOnGrid(const matrix box,
124 const rvec lowerCorner,
125 const rvec upperCorner,
126 const gmx::UpdateGroupsCog* updateGroupsCog,
127 gmx::Range<int> atomRange,
129 gmx::ArrayRef<const int64_t> atomInfo,
130 gmx::ArrayRef<const gmx::RVec> x,
133 nbnxn_atomdata_t* nbat);
135 //! Returns the domain setup
136 DomainSetup domainSetup() const { return domainSetup_; }
138 //! Returns the total number of atoms in the grid set, including padding
139 int numGridAtomsTotal() const { return grids_.back().atomIndexEnd(); }
141 //! Returns the number of local real atoms, i.e. without padded atoms
142 int numRealAtomsLocal() const { return numRealAtomsLocal_; }
144 //! Returns the number of total real atoms, i.e. without padded atoms
145 int numRealAtomsTotal() const { return numRealAtomsTotal_; }
147 //! Returns the atom order on the grid for the local atoms
148 gmx::ArrayRef<const int> getLocalAtomorder() const
150 /* Return the atom order for the home cell (index 0) */
151 const int numIndices = grids_[0].atomIndexEnd() - grids_[0].firstAtomInColumn(0);
153 return gmx::constArrayRefFromArray(atomIndices().data(), numIndices);
156 //! Sets the order of the local atoms to the order grid atom ordering
157 void setLocalAtomOrder();
159 //! Returns the list of grids
160 gmx::ArrayRef<const Grid> grids() const { return grids_; }
162 //! Returns the grid atom indices covering all grids
163 gmx::ArrayRef<const int> cells() const { return gridSetData_.cells; }
165 //! Returns the grid atom indices covering all grids
166 gmx::ArrayRef<const int> atomIndices() const { return gridSetData_.atomIndices; }
168 //! Returns whether we have perturbed non-bonded interactions
169 bool haveFep() const { return haveFep_; }
171 //! Returns the unit cell in \p box
172 void getBox(matrix box) const { copy_mat(box_, box); }
174 //! Returns the maximum number of columns across all grids
175 int numColumnsMax() const { return numColumnsMax_; }
177 //! Sets the maximum number of columns across all grids
178 void setNumColumnsMax(int numColumnsMax) { numColumnsMax_ = numColumnsMax; }
183 DomainSetup domainSetup_;
185 std::vector<Grid> grids_;
186 //! The cell and atom index data which runs over all grids
187 GridSetData gridSetData_;
188 //! Tells whether we have perturbed non-bonded interactions
190 //! The periodic unit-cell
192 //! The number of local real atoms, i.e. without padded atoms, local atoms: 0 to numAtomsLocal_
193 int numRealAtomsLocal_;
194 //! The total number of real atoms, i.e. without padded atoms
195 int numRealAtomsTotal_;
196 //! Working data for constructing a single grid, one entry per thread
197 std::vector<GridWork> gridWork_;
198 //! Maximum number of columns across all grids