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

Detailed Description

A cell list used only to identify nearby atom pairs.

An CellList divides the domain owned by a processor, plus a frame containing ghost particles, into a grid of cells, such that the length of each cell in each direction is greater than a specified cutoff distance. The algorithm works with either generalized or Cartesian coordinates, if used consistently.

All operations of this class are local (no MPI).

Building a 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 local atoms.
AtomStorage::AtomIterator atomIter;
for (storage.begin(atomIter); atomIter.notEnd(); ++atomIter) {
cellList.placeAtom(*atomIter);
}
// Place all ghost atoms
AtomStorage::GhostIterator ghostIter;
for (storage.begin(ghostIter); ghostIter.notEnd(); ++ghostIter){
cellList.placeAtom(*ghostIter);
}
// 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 direciton 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 91 of file ddMd/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, int nCellCut=1)
 Allocate memory for this CellList (generalized coordinates). More...
 
void allocate (int atomCapacity, const Vector &lower, const Vector &upper, double cutoff, int nCellCut=1)
 Allocate memory for this CellList (Cartesian coordinates). More...
 
void makeGrid (const Vector &lower, const Vector &upper, const Vector &cutoffs, int nCellCut=1)
 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

DdMd::CellList::CellList ( )

Constructor.

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

References Util::Dimension.

DdMd::CellList::~CellList ( )
virtual

Destructor.

Deallocates array of Cell objects, if necessary.

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

Member Function Documentation

void DdMd::CellList::allocate ( int  atomCapacity,
const Vector lower,
const Vector upper,
const Vector cutoffs,
int  nCellCut = 1 
)

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
nCellCutnumber of cells per cutoff length

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

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

Referenced by DdMd::PairPotential::save().

void DdMd::CellList::allocate ( int  atomCapacity,
const Vector lower,
const Vector upper,
double  cutoff,
int  nCellCut = 1 
)

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
nCellCutnumber of cells per cutoff length
void DdMd::CellList::makeGrid ( const Vector lower,
const Vector upper,
const Vector cutoffs,
int  nCellCut = 1 
)

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
nCellCutnumber of cells per cutoff length

Definition at line 158 of file ddMd/neighbor/CellList.cpp.

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

Referenced by DdMd::PairPotential::buildCellList().

void DdMd::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 expanded domain for nonbonded ghosts, which extends one cutoff length beyond the domain boundaries the domain boundaries in each direction.

Parameters
atomAtom object to be added.

Definition at line 402 of file ddMd/neighbor/CellList.h.

References DdMd::Atom::position().

Referenced by DdMd::PairPotential::buildCellList().

void DdMd::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 277 of file ddMd/neighbor/CellList.cpp.

References Util::Grid::size().

Referenced by DdMd::PairPotential::buildCellList().

void DdMd::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 310 of file ddMd/neighbor/CellList.cpp.

Referenced by DdMd::PairList::build().

void DdMd::CellList::clear ( )

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

Definition at line 256 of file ddMd/neighbor/CellList.cpp.

References Util::Grid::size().

Referenced by DdMd::PairPotential::buildCellList().

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

Return Grid object by const reference.

Definition at line 427 of file ddMd/neighbor/CellList.h.

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

Dimension of each cell in direction i.

Parameters
iCartesian index i = 0, 1, 2

Definition at line 433 of file ddMd/neighbor/CellList.h.

int DdMd::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 383 of file ddMd/neighbor/CellList.h.

References Util::Dimension.

Referenced by isValid().

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

Return pointer to first local cell in linked list.

Definition at line 448 of file ddMd/neighbor/CellList.h.

Referenced by DdMd::PairList::build(), DdMd::PairPotentialImpl< Interaction >::computeEnergy(), and DdMd::PairPotential::nPair().

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

Return a specified cell by const reference.

Parameters
icell index

Definition at line 439 of file ddMd/neighbor/CellList.h.

int DdMd::CellList::nAtom ( ) const

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

Definition at line 320 of file ddMd/neighbor/CellList.cpp.

Referenced by DdMd::PairPotential::buildCellList().

int DdMd::CellList::nReject ( ) const

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

Definition at line 326 of file ddMd/neighbor/CellList.cpp.

Referenced by DdMd::PairPotential::buildCellList().

int DdMd::CellList::atomCapacity ( ) const

Maximum number of atoms for which space is allocated.

Definition at line 340 of file ddMd/neighbor/CellList.cpp.

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

int DdMd::CellList::cellCapacity ( ) const

Number of cells for which space has been allocated.

Definition at line 346 of file ddMd/neighbor/CellList.cpp.

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

Has memory been allocated for this CellList?

Definition at line 454 of file ddMd/neighbor/CellList.h.

Referenced by isValid().

bool DdMd::CellList::isBuilt ( ) const

Has this CellList been built?

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

Definition at line 352 of file ddMd/neighbor/CellList.cpp.

bool DdMd::CellList::isValid ( ) const

Return true if valid, or throw Exception.

Returns
true if valid, otherwise throw an Exception.

Definition at line 358 of file ddMd/neighbor/CellList.cpp.

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

Referenced by DdMd::PairPotential::buildCellList().


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