Simpatico  v1.10
List of all members | Classes | Public Member Functions
Tools::CellList Class Reference

Detailed Description

A cell list used only to identify nearby atom pairs.

An CellList divides the simulation unit cell into a grid of cells, such that the length of each cell in each direction is greater than a specified cutoff distance.

Building an CellList (Cartesian coordinates);

AtomStorage storage;
CellList cellList;
Vector lower; // Vector of lower bounds (local atoms)
Vector upper; // Vector of upper bounds (local atoms)
Vector cutoffs; // Vector of cutoff lengths for each axis.
int atomCapacity // max number of atoms on this processor
// Set elements of cutoffs vector to same value
for (int i = 0; i < Dimension; ++i) {
cutoffs[i] = cutoff;
}
// Bounds on lower and upper used here to allocate memory.
cellList.allocate(atomCapacity, lower, upper, cutoffs);
// Make the actual grid and clear it.
cellList.makeGrid(lower, upper, cutoffs);
cellList.clear();
// Place all atoms.
AtomStorage::AtomIterator atomIter;
for (storage.begin(atomIter); atomIter.notEnd(); ++atomIter) {
cellList.placeAtom(*atomIter);
}
// Build cell list
cellList.build();

The atomCapacity parameter should be set equal to the sum of atomCapacity and the ghostCapacity of the associated atomStorage, which is the maximum total number of atoms that can exist on this processor.

If the upper and lower bounds and atom coordinates are all expressed in generalized coordinates, which span 0.0 - 1.0 over the primitive periodic cell in each direction, each element of the cutoffs vector is given by a ratio cutoffs[i] = cutoff/length[i], where length[i] is the Cartesian distance across the unit cell along a direction parallel to reciprocal basis vector i.

See Cell documentation for an example of how to iterate over local cells and neighboring atom pairs.

Definition at line 83 of file tools/neighbor/CellList.h.

#include <CellList.h>

Public Member Functions

 CellList ()
 Constructor. More...
 
virtual ~CellList ()
 Destructor. More...
 
void allocate (int atomCapacity, const Vector &lower, const Vector &upper, const Vector &cutoffs)
 Allocate memory for this CellList (generalized coordinates). More...
 
void allocate (int atomCapacity, const Vector &lower, const Vector &upper, double cutoff)
 Allocate memory for this CellList (Cartesian coordinates). More...
 
void makeGrid (const Vector &lower, const Vector &upper, const Vector &cutoffs)
 Make the cell grid (using generalized coordinates). More...
 
void placeAtom (Atom &atom)
 Determine the appropriate cell for an Atom, based on its position. More...
 
void build ()
 Build the cell list. More...
 
void update ()
 Update the cell list. More...
 
void clear ()
 Reset the cell list to its empty state (no atoms). More...
 
const Gridgrid () const
 Return Grid object by const reference. More...
 
double cellLength (int i) const
 Dimension of each cell in direction i. More...
 
int cellIndexFromPosition (const Vector &position) const
 Return the index of the cell that contains a position Vector. More...
 
const Cellbegin () const
 Return pointer to first local cell in linked list. More...
 
const Cellcell (int i) const
 Return a specified cell by const reference. More...
 
int nAtom () const
 Get total number of atoms (local and ghost) in this CellList. More...
 
int nReject () const
 Get number of atoms that were rejected (not placed in cells) More...
 
int atomCapacity () const
 Maximum number of atoms for which space is allocated. More...
 
int cellCapacity () const
 Number of cells for which space has been allocated. More...
 
bool isAllocated () const
 Has memory been allocated for this CellList? More...
 
bool isBuilt () const
 Has this CellList been built? More...
 
bool isValid () const
 Return true if valid, or throw Exception. More...
 

Constructor & Destructor Documentation

Tools::CellList::CellList ( )

Constructor.

Definition at line 21 of file tools/neighbor/CellList.cpp.

References Util::Dimension.

Tools::CellList::~CellList ( )
virtual

Destructor.

Definition at line 38 of file tools/neighbor/CellList.cpp.

Member Function Documentation

void Tools::CellList::allocate ( int  atomCapacity,
const Vector lower,
const Vector upper,
const Vector cutoffs 
)

Allocate memory for this CellList (generalized coordinates).

This function:

  • Allocates an array of atomCapacity CellList::Tag objects.
  • Allocates an array of atomCapacity CellAtom objects.
  • Allocates an array of Cell objects sized for this boundary.

The elements of the lower, upper, and cutoffs parameters should contain the lower and upper coordinate bounds for this processor, and cutoff values in each direction, for a boundary was chosen to be larger than any that will be encountered during the simulation These parameters are used only to allocate memory.

This version of the function is designed for use with generalized coordinates. See the makeGrid() method for a discussion of the parameter values in generalized coordinates.

Parameters
atomCapacitydimension of global array of atoms
lowerlower coordinate bounds for this processor
upperupper coordinate bounds for this processor
cutoffspair cutoff distance in each direction

Definition at line 44 of file tools/neighbor/CellList.cpp.

References Util::DArray< Data >::allocate(), Util::Dimension, Util::Grid::dimension(), Util::Grid::rank(), Util::Grid::setDimensions(), Tools::Cell::setLastCell(), Tools::Cell::setNextCell(), Util::Grid::size(), and UTIL_THROW.

Referenced by Tools::PairEnergy::setup().

void Tools::CellList::allocate ( int  atomCapacity,
const Vector lower,
const Vector upper,
double  cutoff 
)

Allocate memory for this CellList (Cartesian coordinates).

This function is designed for use with Cartesian coordinates, for which the lower, upper and cutoff parameters all have dimensions of length. The function calls the allocate() method with a Vector of cutoffs internally, after setting every element of the Vector to the same value.

Parameters
atomCapacitydimension of global array of atoms
lowerlower bound for this processor in maximum boundary
upperupper bound for this processor in maximum boundary
cutoffpair cutoff distance in each direction
void Tools::CellList::makeGrid ( const Vector lower,
const Vector upper,
const Vector cutoffs 
)

Make the cell grid (using generalized coordinates).

This method makes a cell grid in which the number of cells in each direction i is chosen such that the dimension of each cell that direction is greater than or equal to cutoff[i].

The elements of lower and upper should be upper and lower bounds for coordinates of local atoms on this processor, in generalized coordinates in which the entire periodic unit cell spans 0.0 to 1.0. On input, each element of cutoff[i] should be equal to the minimum cell length in direction i in generalized coordinates. This is given by the ratio cutoff[i] = pairCutoff/length[i], where pairCutoff is the maximum range of nonbonded interactions, and length[i] is the distance across the primitive unit cell along the direction parallel to reciprocal lattice basis vector i.

Parameters
lowerlower bound of local atom coordinates.
upperupper bound of local atom coordinates.
cutoffspair cutoff length in each direction

Definition at line 146 of file tools/neighbor/CellList.cpp.

References Util::FSArray< Data, Capacity >::append(), Util::Dimension, and Util::Grid::size().

Referenced by Tools::PairEnergy::setup().

void Tools::CellList::placeAtom ( Atom atom)
inline

Determine the appropriate cell for an Atom, based on its position.

This method does not place the atom in a cell, but calculates a cell index and retains the value, which is used to place atoms in the build() method.

The method quietly does nothing if the atom is outside the domain.

Parameters
atomAtom object to be added.

Definition at line 413 of file tools/neighbor/CellList.h.

References Tools::Atom::position.

Referenced by Tools::PairEnergy::sample().

void Tools::CellList::build ( )

Build the cell list.

This method must be called after completing a loop over atoms, in which placeAtom() is called for each atom. The build() method uses information gathered in this loop to build and fill all of the cells.

Definition at line 192 of file tools/neighbor/CellList.cpp.

References Util::Grid::size().

Referenced by Tools::PairEnergy::sample().

void Tools::CellList::update ( )

Update the cell list.

This method must be called after the cell list is built, and after transformation back to Cartesian coordinates.

Definition at line 225 of file tools/neighbor/CellList.cpp.

Referenced by Tools::PairEnergy::sample().

void Tools::CellList::clear ( )

Reset the cell list to its empty state (no atoms).

Definition at line 171 of file tools/neighbor/CellList.cpp.

References Util::Grid::size().

Referenced by Tools::PairEnergy::sample().

const Grid & Tools::CellList::grid ( ) const
inline

Return Grid object by const reference.

Definition at line 438 of file tools/neighbor/CellList.h.

double Tools::CellList::cellLength ( int  i) const
inline

Dimension of each cell in direction i.

Parameters
iCartesian index i = 0, 1, 2

Definition at line 444 of file tools/neighbor/CellList.h.

int Tools::CellList::cellIndexFromPosition ( const Vector position) const
inline

Return the index of the cell that contains a position Vector.

Returns a null value of -1 if the position is outside the expanded domain for nonbonded ghost atoms.

Parameters
positionposition Vector, inside the boundary

Definition at line 394 of file tools/neighbor/CellList.h.

References Util::Dimension.

Referenced by isValid().

const Cell * Tools::CellList::begin ( ) const
inline

Return pointer to first local cell in linked list.

Definition at line 456 of file tools/neighbor/CellList.h.

Referenced by Tools::PairEnergy::sample().

const Cell & Tools::CellList::cell ( int  i) const
inline

Return a specified cell by const reference.

Parameters
icell index

Definition at line 450 of file tools/neighbor/CellList.h.

int Tools::CellList::nAtom ( ) const

Get total number of atoms (local and ghost) in this CellList.

Definition at line 235 of file tools/neighbor/CellList.cpp.

int Tools::CellList::nReject ( ) const

Get number of atoms that were rejected (not placed in cells)

Definition at line 241 of file tools/neighbor/CellList.cpp.

int Tools::CellList::atomCapacity ( ) const

Maximum number of atoms for which space is allocated.

Definition at line 255 of file tools/neighbor/CellList.cpp.

References Util::Array< Data >::capacity().

int Tools::CellList::cellCapacity ( ) const

Number of cells for which space has been allocated.

Definition at line 261 of file tools/neighbor/CellList.cpp.

bool Tools::CellList::isAllocated ( ) const
inline

Has memory been allocated for this CellList?

Definition at line 462 of file tools/neighbor/CellList.h.

Referenced by isValid().

bool Tools::CellList::isBuilt ( ) const

Has this CellList been built?

IsBuilt is set true by the build() method and false by clear().

Definition at line 267 of file tools/neighbor/CellList.cpp.

bool Tools::CellList::isValid ( ) const

Return true if valid, or throw Exception.

Returns
true if valid, otherwise throw an Exception.

Definition at line 273 of file tools/neighbor/CellList.cpp.

References Tools::Cell::atomCapacity(), Tools::Cell::atomPtr(), Util::Array< Data >::capacity(), cellIndexFromPosition(), isAllocated(), Tools::Cell::nAtom(), Util::Grid::size(), and UTIL_THROW.

Referenced by Tools::PairEnergy::sample().


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