Gromacs  2026.0-dev-20250116-fa3fd9d
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
List of all members | Classes | Public Member Functions | Static Public Member Functions
gmx::Grid Class Reference

#include <gromacs/nbnxm/grid.h>

Description

A pair-search grid object for one domain decomposition zone.

This is a rectangular 3D grid covering a potentially non-rectangular volume which is either the whole unit cell or the local zone or part of a non-local zone when using domain decomposition. The grid cells are even spaced along x/y and irregular along z. Each cell is sub-divided into atom clusters. With a CPU geometry, each cell contains 1 or 2 clusters. With a GPU geometry, each cell contains up to 8 clusters. The geometry is set by the pairlist type which is the only argument of the constructor.

When multiple grids are used, i.e. with domain decomposition, we want to avoid the overhead of multiple coordinate arrays or extra indexing. Therefore each grid stores a cell offset, so a contiguous cell index can be used to index atom arrays. All methods returning atom indices return indices which index into a full atom array.

Note that when atom groups, instead of individual atoms, are assigned to grid cells, individual atoms can be geometrically outside the cell and grid that they have been assigned to (as determined by the center or geometry of the atom group they belong to).

Classes

struct  Geometry
 The cluster and cell geometry of a grid. More...
 

Public Member Functions

 Grid (PairlistType pairlistType, int ddZone, const bool &haveFep, PinningPolicy pinningPolicy)
 Constructs a grid given the type of pairlist.
 
const Geometrygeometry () const
 Returns the geometry of the grid cells.
 
int ddZone () const
 Returns the zone this grid belong to.
 
const GridDimensionsdimensions () const
 Returns the dimensions of the grid.
 
int numColumns () const
 Returns the total number of grid columns.
 
int numCells () const
 Returns the total number of grid cells.
 
int cellOffset () const
 Returns the cell offset of (the first cell of) this grid in the list of cells combined over all grids.
 
int numCellsColumnMax () const
 Returns the maximum number of grid cells in a column.
 
int srcAtomBegin () const
 Returns the start of the source atom range mapped to this grid.
 
int srcAtomEnd () const
 Returns the end of the source atom range mapped to this grid.
 
int firstCellInColumn (int columnIndex) const
 Returns the first cell index in the grid, starting at 0 in this grid.
 
int numCellsInColumn (int columnIndex) const
 Returns the number of cells in the column.
 
int firstAtomInColumn (int columnIndex) const
 Returns the index of the first atom in the column.
 
int numAtomsInColumn (int columnIndex) const
 Returns the number of real atoms in the column.
 
ArrayRef< const int > cxy_na () const
 Returns a view of the number of non-filler, atoms for each grid column. More...
 
ArrayRef< const int > cxy_ind () const
 Returns a view of the grid-local cell index for each grid column. More...
 
int numAtomsPerCell () const
 Returns the number of real atoms in the column.
 
int paddedNumAtomsInColumn (int columnIndex) const
 Returns the number of atoms in the column including padding.
 
int atomIndexEnd () const
 Returns the end of the atom index range on the grid, including padding.
 
bool clusterIsPerturbed (int clusterIndex) const
 Returns whether any atom in the cluster is perturbed.
 
bool atomIsPerturbed (int clusterIndex, int atomIndexInCluster) const
 Returns whether the given atom in the cluster is perturbed.
 
unsigned int fepBits (int clusterIndex) const
 Returns the free-energy perturbation bits for the cluster.
 
ArrayRef< const BoundingBoxiBoundingBoxes () const
 Returns the i-bounding boxes for all clusters on the grid.
 
ArrayRef< const BoundingBoxjBoundingBoxes () const
 Returns the j-bounding boxes for all clusters on the grid.
 
ArrayRef< const float > packedBoundingBoxes () const
 Returns the packed bounding boxes for all clusters on the grid, empty with a CPU list.
 
ArrayRef< const BoundingBox1DzBoundingBoxes () const
 Returns the bounding boxes along z for all i-cells on the grid.
 
ArrayRef< const int > clusterFlags () const
 Returns the flags for all clusters on the grid.
 
ArrayRef< const int > numClustersPerCell () const
 Returns the number of clusters for all cells on the grid, empty with a CPU geometry.
 
int atomToCluster (int atomIndex) const
 Returns the cluster index for an atom.
 
int numClusters () const
 Returns the total number of clusters on the grid.
 
void resizeBoundingBoxesAndFlags (const int maxNumCells)
 Resizes the bouding box and FEP flag lists for at most maxNumCells.
 
void setDimensions (int ddZone, int numAtomsTotal, int numAtomsWithoutFillers, const RVec &lowerCorner, const RVec &upperCorner, real *atomDensity, real maxAtomGroupRadius)
 Sets the grid dimensions. More...
 
void setCellIndices (int ddZone, int cellOffset, GridSetData *gridSetData, ArrayRef< GridWork > gridWork, Range< int > atomRange, ArrayRef< const int32_t > atomInfo, ArrayRef< const RVec > x, nbnxn_atomdata_t *nbat)
 Sets the cell indices using indices in gridSetData and gridWork.
 
void setNonLocalGrid (int ddZone, const GridDimensions &dimensions, ArrayRef< const std::pair< int, int >> columns, int cellOffset, ArrayRef< const int32_t > atomInfo, ArrayRef< const RVec > x, GridSetData *gridSetData, nbnxn_atomdata_t *nbat)
 Sets a non-local grid using data communicated from a different domain. More...
 

Static Public Member Functions

static void calcColumnIndices (const GridDimensions &gridDims, const UpdateGroupsCog *updateGroupsCog, Range< int > atomRange, ArrayRef< const RVec > x, int dd_zone, const int *move, int thread, int nthread, ArrayRef< int > cell, ArrayRef< int > cxy_na)
 Determine in which grid columns atoms should go, store cells and atom counts in cell and cxy_na.
 

Member Function Documentation

ArrayRef<const int> gmx::Grid::cxy_ind ( ) const
inline

Returns a view of the grid-local cell index for each grid column.

Todo:
Needs a useful name.
ArrayRef<const int> gmx::Grid::cxy_na ( ) const
inline

Returns a view of the number of non-filler, atoms for each grid column.

Todo:
Needs a useful name.
void gmx::Grid::setDimensions ( int  ddZone,
int  numAtomsTotal,
int  numAtomsWithoutFillers,
const RVec lowerCorner,
const RVec upperCorner,
real atomDensity,
real  maxAtomGroupRadius 
)

Sets the grid dimensions.

Parameters
[in]ddZoneThe domain decomposition zone index
[in]numAtomsTotalThe total number of atoms to put onto this grid
[in]numAtomsWithoutFillersThe number of atoms that are not filler particles
[in]lowerCornerThe minimum Cartesian coordinates of the grid
[in]upperCornerThe maximum Cartesian coordinates of the grid
[in,out]atomDensityThe atom density, will be computed when <= 0
[in]maxAtomGroupRadiusThe maximum radius of atom groups
void gmx::Grid::setNonLocalGrid ( int  ddZone,
const GridDimensions dimensions,
ArrayRef< const std::pair< int, int >>  columns,
int  cellOffset,
ArrayRef< const int32_t >  atomInfo,
ArrayRef< const RVec x,
GridSetData gridSetData,
nbnxn_atomdata_t nbat 
)

Sets a non-local grid using data communicated from a different domain.

Parameters
[in]ddZoneThe domain decomposition zone this grid belongs to
[in]dimensionsThe dimensions of the grid
[in]columnsA list of column indices and number of cells for a column, the list should be ordered on column index
[in]cellOffsetThe offset of this grid in the list of cells over all grids
[in]atomInfoA list of information for all local and non-local atoms
[in]xThe coordinates for all local and non-local atoms
[in,out]gridSetDataThe data shared over all grids
[in,out]nbatThe NBNxM atom data, used here for storing the atom coordinates

The documentation for this class was generated from the following files: