2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2019,2020, 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 Grid class.
41 * This class provides functionality for setting up and accessing atoms
42 * on a grid for one domain decomposition zone. This grid is used for
43 * generating cluster pair lists for computing non-bonded pair interactions.
44 * The grid consists of a regular array of columns along dimensions x and y.
45 * Along z the number of cells and their boundaries vary between the columns.
46 * Each cell can hold one or more clusters of atoms, depending on the grid
47 * geometry, which is set by the pair-list type.
49 * \author Berk Hess <hess@kth.se>
50 * \ingroup module_nbnxm
53 #ifndef GMX_NBNXM_GRID_H
54 #define GMX_NBNXM_GRID_H
59 #include "gromacs/gpu_utils/hostallocator.h"
60 #include "gromacs/math/vectypes.h"
61 #include "gromacs/utility/alignedallocator.h"
62 #include "gromacs/utility/arrayref.h"
63 #include "gromacs/utility/range.h"
65 struct gmx_domdec_zones_t
;
66 struct nbnxn_atomdata_t
;
68 enum class PairlistType
;
72 class UpdateGroupsCog
;
82 * \brief Bounding box for a nbnxm atom cluster
84 * \note Should be aligned in memory to enable 4-wide SIMD operations.
89 * \brief Corner for the bounding box, padded with one element to enable 4-wide SIMD operations
93 //! Returns a corner with the minimum coordinates along each dimension
94 static Corner
min(const Corner
& c1
, const Corner
& c2
)
98 cMin
.x
= std::min(c1
.x
, c2
.x
);
99 cMin
.y
= std::min(c1
.y
, c2
.y
);
100 cMin
.z
= std::min(c1
.z
, c2
.z
);
101 /* This value of the padding is irrelevant, as long as it
102 * is initialized. We use min to allow auto-vectorization.
104 cMin
.padding
= std::min(c1
.padding
, c2
.padding
);
109 //! Returns a corner with the maximum coordinates along each dimension
110 static Corner
max(const Corner
& c1
, const Corner
& c2
)
114 cMax
.x
= std::max(c1
.x
, c2
.x
);
115 cMax
.y
= std::max(c1
.y
, c2
.y
);
116 cMax
.z
= std::max(c1
.z
, c2
.z
);
117 cMax
.padding
= std::max(c1
.padding
, c2
.padding
);
122 //! Returns a pointer for SIMD loading of a Corner object
123 const float* ptr() const { return &x
; }
125 //! Returns a pointer for SIMD storing of a Corner object
126 float* ptr() { return &x
; }
134 //! padding, unused, but should be set to avoid operations on unitialized data
138 //! lower, along x and y and z, corner
140 //! upper, along x and y and z, corner
145 * \brief Bounding box for one dimension of a grid cell
161 * \brief A pair-search grid object for one domain decomposition zone
163 * This is a rectangular 3D grid covering a potentially non-rectangular
164 * volume which is either the whole unit cell or the local zone or part
165 * of a non-local zone when using domain decomposition. The grid cells
166 * are even spaced along x/y and irregular along z. Each cell is sub-divided
167 * into atom clusters. With a CPU geometry, each cell contains 1 or 2 clusters.
168 * With a GPU geometry, each cell contains up to 8 clusters. The geometry is
169 * set by the pairlist type which is the only argument of the constructor.
171 * When multiple grids are used, i.e. with domain decomposition, we want
172 * to avoid the overhead of multiple coordinate arrays or extra indexing.
173 * Therefore each grid stores a cell offset, so a contiguous cell index
174 * can be used to index atom arrays. All methods returning atom indices
175 * return indices which index into a full atom array.
177 * Note that when atom groups, instead of individual atoms, are assigned
178 * to grid cells, individual atoms can be geometrically outside the cell
179 * and grid that they have been assigned to (as determined by the center
180 * or geometry of the atom group they belong to).
186 * \brief The cluster and cell geometry of a grid
190 //! Constructs the cluster/cell geometry given the type of pairlist
191 Geometry(PairlistType pairlistType
);
193 //! Is this grid simple (CPU) or hierarchical (GPU)
195 //! Number of atoms per cluster
196 int numAtomsICluster
;
197 //! Number of atoms for list j-clusters
198 int numAtomsJCluster
;
199 //! Number of atoms per cell
202 int numAtomsICluster2Log
;
205 //! The physical dimensions of a grid \internal
208 //! The lower corner of the (local) grid
210 //! The upper corner of the (local) grid
212 //! The physical grid size: upperCorner - lowerCorner
214 //! An estimate for the atom number density of the region targeted by the grid
216 //! The maximum distance an atom can be outside of a cell and outside of the grid
217 real maxAtomGroupRadius
;
218 //! Size of cell along dimension x and y
219 real cellSize
[DIM
- 1];
220 //! 1/size of a cell along dimensions x and y
221 real invCellSize
[DIM
- 1];
222 //! The number of grid cells along dimensions x and y
223 int numCells
[DIM
- 1];
226 //! Constructs a grid given the type of pairlist
227 Grid(PairlistType pairlistType
, const bool& haveFep
);
229 //! Returns the geometry of the grid cells
230 const Geometry
& geometry() const { return geometry_
; }
232 //! Returns the dimensions of the grid
233 const Dimensions
& dimensions() const { return dimensions_
; }
235 //! Returns the total number of grid columns
236 int numColumns() const { return dimensions_
.numCells
[XX
] * dimensions_
.numCells
[YY
]; }
238 //! Returns the total number of grid cells
239 int numCells() const { return numCellsTotal_
; }
241 //! Returns the cell offset of (the first cell of) this grid in the list of cells combined over all grids
242 int cellOffset() const { return cellOffset_
; }
244 //! Returns the maximum number of grid cells in a column
245 int numCellsColumnMax() const { return numCellsColumnMax_
; }
247 //! Returns the start of the source atom range mapped to this grid
248 int srcAtomBegin() const { return srcAtomBegin_
; }
250 //! Returns the end of the source atom range mapped to this grid
251 int srcAtomEnd() const { return srcAtomEnd_
; }
253 //! Returns the first cell index in the grid, starting at 0 in this grid
254 int firstCellInColumn(int columnIndex
) const { return cxy_ind_
[columnIndex
]; }
256 //! Returns the number of cells in the column
257 int numCellsInColumn(int columnIndex
) const
259 return cxy_ind_
[columnIndex
+ 1LL] - cxy_ind_
[columnIndex
];
262 //! Returns the index of the first atom in the column
263 int firstAtomInColumn(int columnIndex
) const
265 return (cellOffset_
+ cxy_ind_
[columnIndex
]) * geometry_
.numAtomsPerCell
;
268 //! Returns the number of real atoms in the column
269 int numAtomsInColumn(int columnIndex
) const { return cxy_na_
[columnIndex
]; }
271 /*! \brief Returns a view of the number of non-filler, atoms for each grid column
273 * \todo Needs a useful name. */
274 gmx::ArrayRef
<const int> cxy_na() const { return cxy_na_
; }
275 /*! \brief Returns a view of the grid-local cell index for each grid column
277 * \todo Needs a useful name. */
278 gmx::ArrayRef
<const int> cxy_ind() const { return cxy_ind_
; }
280 //! Returns the number of real atoms in the column
281 int numAtomsPerCell() const { return geometry_
.numAtomsPerCell
; }
283 //! Returns the number of atoms in the column including padding
284 int paddedNumAtomsInColumn(int columnIndex
) const
286 return numCellsInColumn(columnIndex
) * geometry_
.numAtomsPerCell
;
289 //! Returns the end of the atom index range on the grid, including padding
290 int atomIndexEnd() const { return (cellOffset_
+ numCellsTotal_
) * geometry_
.numAtomsPerCell
; }
292 //! Returns whether any atom in the cluster is perturbed
293 bool clusterIsPerturbed(int clusterIndex
) const { return fep_
[clusterIndex
] != 0U; }
295 //! Returns whether the given atom in the cluster is perturbed
296 bool atomIsPerturbed(int clusterIndex
, int atomIndexInCluster
) const
298 return (fep_
[clusterIndex
] & (1 << atomIndexInCluster
)) != 0U;
301 //! Returns the free-energy perturbation bits for the cluster
302 unsigned int fepBits(int clusterIndex
) const { return fep_
[clusterIndex
]; }
304 //! Returns the i-bounding boxes for all clusters on the grid
305 gmx::ArrayRef
<const BoundingBox
> iBoundingBoxes() const { return bb_
; }
307 //! Returns the j-bounding boxes for all clusters on the grid
308 gmx::ArrayRef
<const BoundingBox
> jBoundingBoxes() const { return bbj_
; }
310 //! Returns the packed bounding boxes for all clusters on the grid, empty with a CPU list
311 gmx::ArrayRef
<const float> packedBoundingBoxes() const { return pbb_
; }
313 //! Returns the bounding boxes along z for all cells on the grid
314 gmx::ArrayRef
<const BoundingBox1D
> zBoundingBoxes() const { return bbcz_
; }
316 //! Returns the flags for all clusters on the grid
317 gmx::ArrayRef
<const int> clusterFlags() const { return flags_
; }
319 //! Returns the number of clusters for all cells on the grid, empty with a CPU geometry
320 gmx::ArrayRef
<const int> numClustersPerCell() const { return numClusters_
; }
322 //! Returns the cluster index for an atom
323 int atomToCluster(int atomIndex
) const { return (atomIndex
>> geometry_
.numAtomsICluster2Log
); }
325 //! Returns the total number of clusters on the grid
326 int numClusters() const
328 if (geometry_
.isSimple
)
330 return numCellsTotal_
;
334 return numClustersTotal_
;
338 //! Sets the grid dimensions
339 void setDimensions(int ddZone
,
341 gmx::RVec lowerCorner
,
342 gmx::RVec upperCorner
,
344 real maxAtomGroupRadius
,
346 gmx::PinningPolicy pinningPolicy
);
348 //! Sets the cell indices using indices in \p gridSetData and \p gridWork
349 void setCellIndices(int ddZone
,
351 GridSetData
* gridSetData
,
352 gmx::ArrayRef
<GridWork
> gridWork
,
353 gmx::Range
<int> atomRange
,
355 gmx::ArrayRef
<const gmx::RVec
> x
,
357 nbnxn_atomdata_t
* nbat
);
359 //! Determine in which grid columns atoms should go, store cells and atom counts in \p cell and \p cxy_na
360 static void calcColumnIndices(const Grid::Dimensions
& gridDims
,
361 const gmx::UpdateGroupsCog
* updateGroupsCog
,
362 gmx::Range
<int> atomRange
,
363 gmx::ArrayRef
<const gmx::RVec
> x
,
368 gmx::ArrayRef
<int> cell
,
369 gmx::ArrayRef
<int> cxy_na
);
372 /*! \brief Fill a pair search cell with atoms
374 * Potentially sorts atoms and sets the interaction flags.
376 void fillCell(GridSetData
* gridSetData
,
377 nbnxn_atomdata_t
* nbat
,
381 gmx::ArrayRef
<const gmx::RVec
> x
,
382 BoundingBox gmx_unused
* bb_work_aligned
);
384 //! Spatially sort the atoms within the given column range, for CPU geometry
385 void sortColumnsCpuGeometry(GridSetData
* gridSetData
,
388 gmx::ArrayRef
<const gmx::RVec
> x
,
389 nbnxn_atomdata_t
* nbat
,
390 gmx::Range
<int> columnRange
,
391 gmx::ArrayRef
<int> sort_work
);
393 //! Spatially sort the atoms within the given column range, for GPU geometry
394 void sortColumnsGpuGeometry(GridSetData
* gridSetData
,
397 gmx::ArrayRef
<const gmx::RVec
> x
,
398 nbnxn_atomdata_t
* nbat
,
399 gmx::Range
<int> columnRange
,
400 gmx::ArrayRef
<int> sort_work
);
403 //! The geometry of the grid clusters and cells
405 //! The physical dimensions of the grid
406 Dimensions dimensions_
;
408 //! The total number of cells in this grid
410 //! Index in nbs->cell corresponding to cell 0
412 //! The maximum number of cells in a column
413 int numCellsColumnMax_
;
415 //! The start of the source atom range mapped to this grid
417 //! The end of the source atom range mapped to this grid
421 /*! \brief The number of, non-filler, atoms for each grid column.
423 * \todo Needs a useful name. */
424 gmx::HostVector
<int> cxy_na_
;
425 /*! \brief The grid-local cell index for each grid column
427 * \todo Needs a useful name. */
428 gmx::HostVector
<int> cxy_ind_
;
430 //! The number of cluster for each cell
431 std::vector
<int> numClusters_
;
434 //! Bounding boxes in z for the cells
435 std::vector
<BoundingBox1D
> bbcz_
;
436 //! 3D bounding boxes for the sub cells
437 std::vector
<BoundingBox
, gmx::AlignedAllocator
<BoundingBox
>> bb_
;
438 //! 3D j-bounding boxes for the case where the i- and j-cluster sizes are different
439 std::vector
<BoundingBox
, gmx::AlignedAllocator
<BoundingBox
>> bbjStorage_
;
440 //! 3D j-bounding boxes
441 gmx::ArrayRef
<BoundingBox
> bbj_
;
442 //! 3D bounding boxes in packed xxxx format per cell
443 std::vector
<float, gmx::AlignedAllocator
<float>> pbb_
;
445 //! Tells whether we have perturbed interactions, authorative source is in GridSet (never modified)
446 const bool& haveFep_
;
448 /* Bit-flag information */
449 //! Flags for properties of clusters in each cell
450 std::vector
<int> flags_
;
451 //! Signal bits for atoms in each cell that tell whether an atom is perturbed
452 std::vector
<unsigned int> fep_
;
455 //! Total number of clusters, used for printing
456 int numClustersTotal_
;