Simpatico  v1.10
List of all members | Public Types | Public Member Functions | Static Public Attributes | Friends
McMd::CellList Class Reference

Detailed Description

A cell list for Atom objects in a periodic system boundary.

A CellList divides the volume within a Boundary into a grid of cells (i.e., regions), such that the length of each cell in each of 3 directions is greater * than a specified cutoff distance. Each cell has an integer index between 0 and totCells_ - 1. A CellList can be used to identify all atoms that are within a cutoff distance of any point within a system, as required to calculate nonbonded energies and forces.

In this implementation, each cell is represented by a Cell object. Each cell contains an small array of pointers to all all atoms in the associated volume. This implementation is designed to allow fast addition and removal of individual atoms. This implementation also wastes some memory, compared to a linked list implementation, because the dimension Cell::MaxAtomCell of the array of pointers in each Cell must be large enough to accomodate the maximum number of atoms that will ever be encountered in an individual cell.

CellList is a non-polymorphic class, with no virtual functions, and a non-virtual destructor. Do not derive subclasses from it.

Definition at line 56 of file mcMd/neighbor/CellList.h.

#include <CellList.h>

Public Types

typedef FSArray< Atom *, MaxNeighborNeighborArray
 Static array for holding neighbors in a cell list. More...
 

Public Member Functions

 CellList ()
 Constructor. More...
 
virtual ~CellList ()
 Destructor. More...
 
void setAtomCapacity (int atomCapacity)
 Set atom capacity. More...
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 Serialize to/from an Archive. More...
 
void setup (const Boundary &boundary, double cutoff)
 Setup grid of empty cells. More...
 
void clear ()
 Sets all Cell objects to empty state (no Atoms). More...
 
int cellIndexFromPosition (const Vector &pos) const
 Return index of the cell that contains the position pos = (x,y,z). More...
 
void addAtom (Atom &atom)
 Add a Atom to the appropriate cell, based on its position. More...
 
void deleteAtom (Atom &atom)
 Delete a Atom object from its cell. More...
 
void updateAtomCell (Atom &atom, const Vector &pos)
 Update CellList to reflect a new atom position. More...
 
void getNeighbors (const Vector &pos, NeighborArray &neighbors) const
 Fill a NeighborArray with pointers to atoms near a specified position. More...
 
void getCellNeighbors (int ic, NeighborArray &neighbors, int &nInCell) const
 Fill an array with pointers to atoms in a cell and neighboring cells. More...
 
int gridDimension (int i) const
 Number of cells along axis i. More...
 
int totCells () const
 Get total number of cells in this CellList. More...
 
int nAtom () const
 Get total number of atoms in this CellList. More...
 
int atomCapacity () const
 Get the maximum allowed atom index + 1. More...
 
bool isValid (int nAtom=-1) const
 Return true if valid, or throw Exception. More...
 

Static Public Attributes

static const int MaxNeighbor = 27*Cell::MaxAtomCell
 Maximum possible number of neighboring atoms. More...
 

Friends

class ::CellListTest
 Grant access to unit test class. More...
 

Constructor & Destructor Documentation

McMd::CellList::CellList ( )

Constructor.

Definition at line 20 of file mcMd/neighbor/CellList.cpp.

References Util::Dimension.

McMd::CellList::~CellList ( )
virtual

Destructor.

Deallocates array of Cell objects, if necessary.

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

Member Function Documentation

void McMd::CellList::setAtomCapacity ( int  atomCapacity)

Set atom capacity.

This function sets the atom capacity (the maximum number of atoms), and allocates an array of atomCapacity CellTag objects.

Parameters
atomCapacitydimension of global array of atoms

Definition at line 93 of file mcMd/neighbor/CellList.cpp.

References UTIL_THROW.

Referenced by McMd::ClusterIdentifier::initialize(), McMd::PairList::initialize(), McMd::McPairPotentialImpl< Interaction >::loadParameters(), McMd::McPairPotentialImpl< Interaction >::readParameters(), McMd::Crosslinker::setup(), and McMd::Generator::setupCellList().

template<class Archive >
void McMd::CellList::serialize ( Archive &  ar,
const unsigned int  version 
)

Serialize to/from an Archive.

Parameters
ararchive
versionarchive version id

Definition at line 473 of file mcMd/neighbor/CellList.h.

References Util::DArray< Data >::allocate().

void McMd::CellList::setup ( const Boundary boundary,
double  cutoff 
)

Setup grid of empty cells.

The number of cells in each direction is chosen such that the dimension of each cell in each direction is greater than or equal to the cutoff parameter. To calculate nonbonded pair interaction energies, the cutoff parameter should thus be equal to or greater than the maximum range of nonbonded interactions.

After calculating the grid dimensions, this this function may allocate or resizes the array of cells if necessary, and then calls clear(). It also retains a pointer to the Boundary, which is used to determine the correct cell by the cellIndexFromPosition() function.

Parameters
boundaryBoundary object for the system
cutoffMinimum dimension of a cell in any direction

Definition at line 114 of file mcMd/neighbor/CellList.cpp.

References Simp::OrthorhombicBoundary::lengths(), and UTIL_THROW.

Referenced by McMd::McPairPotential::buildCellList(), McMd::ClusterIdentifier::identifyClusters(), McMd::Crosslinker::setup(), McMd::PairList::setup(), and McMd::Generator::setupCellList().

void McMd::CellList::clear ( )

Sets all Cell objects to empty state (no Atoms).

This function does not change the cell grid dimensions, and does not allocate or reallocate any memory.

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

Referenced by McMd::PairList::clear(), and McMd::Crosslinker::sample().

int McMd::CellList::cellIndexFromPosition ( const Vector pos) const
inline

Return index of the cell that contains the position pos = (x,y,z).

Parameters
posposition vector, as an array pos = {x, y, z}
Returns
index of cell that contains position pos

Definition at line 361 of file mcMd/neighbor/CellList.h.

References Simp::OrthorhombicBoundary::transformCartToGen().

void McMd::CellList::addAtom ( Atom atom)
inline

Add a Atom to the appropriate cell, based on its position.

Parameters
atomAtom object to be added.

Definition at line 397 of file mcMd/neighbor/CellList.h.

References McMd::Atom::id(), and McMd::Atom::position().

Referenced by McMd::Generator::attemptPlaceAtom(), McMd::McPairPotential::buildCellList(), McMd::ClusterIdentifier::identifyClusters(), and McMd::Crosslinker::sample().

void McMd::CellList::deleteAtom ( Atom atom)
inline

Delete a Atom object from its cell.

Parameters
atomAtom object to be deleted

Definition at line 386 of file mcMd/neighbor/CellList.h.

References McMd::CellTag::cellId, and McMd::Atom::id().

Referenced by McMd::LinearGenerator::attemptPlaceMolecule().

void McMd::CellList::updateAtomCell ( Atom atom,
const Vector pos 
)
inline

Update CellList to reflect a new atom position.

If a new atom position pos lies in a new cell, with a cell index icell != atom.cellId, then delete the atom from the old Cell object and add it to the new one. If the new cell is the same as the old one (icell == atom.CellId), do nothing and return.

This function updates the CellList, but does not update the coordinates stored in the pos member of the associated Atom object. To update both, use the System::moveAtom() function.

Parameters
atomAtom object to be moved
posnew Atom position

Definition at line 409 of file mcMd/neighbor/CellList.h.

References McMd::CellTag::cellId, and McMd::Atom::id().

void McMd::CellList::getNeighbors ( const Vector pos,
NeighborArray neighbors 
) const
inline

Fill a NeighborArray with pointers to atoms near a specified position.

Upon return, array neighbors contains pointers to all Atom objects in the cell containing the position pos, and in all neighboring cells. The number of such atoms is given by the logical size() of this array.

See the definition of the typedef CellList::NeighborArray.

Parameters
posposition Vector, pos = {x, y, z}
neighborsarray of pointers to neighboring Atoms

Definition at line 426 of file mcMd/neighbor/CellList.h.

Referenced by McMd::McPairPotentialImpl< Interaction >::atomEnergy(), McMd::Generator::attemptPlaceAtom(), McMd::ClusterIdentifier::initialize(), McMd::McPairPotentialImpl< Interaction >::moleculeEnergy(), McMd::SliplinkMove::move(), McMd::SliplinkerAll::move(), McMd::GcSliplinkMove::move(), McMd::SliplinkerEnd::move(), McMd::McMuExchange::sample(), and McMd::RingOctaRebridgeMove::scanBridge().

void McMd::CellList::getCellNeighbors ( int  ic,
NeighborArray neighbors,
int &  nInCell 
) const

Fill an array with pointers to atoms in a cell and neighboring cells.

Upon return, the NeighborArray neighbors contains pointers to all of the atoms in cell number ic and all neighboring cells. The number of such Atoms is given by the size() of this FSArray < Atom* > array. Also, upon return, nInCell is the number of Atoms in cell ic. Pointers to the atoms in cell ic are listed in the elements neighbors[0], ... , neighbors[nInCell-1] of the neighbors array, while those in neighboring are in elements neighbors[nInCell], ... , neighbors[neighbors.size()-1].

Parameters
iccell index
neighborsFSArray containing pointers to neighbor Atoms
nInCellnumber of Atoms in cell ic.

Definition at line 147 of file mcMd/neighbor/CellList.cpp.

References Util::FSArray< Data, Capacity >::append(), McMd::Cell::atomPtr(), Util::FSArray< Data, Capacity >::clear(), and McMd::Cell::firstClearPos().

Referenced by McMd::PairList::build(), McMd::McPairPotentialImpl< Interaction >::computeEnergy(), McMd::McPairPerturbation< Interaction >::derivative(), McMd::McPairExternalPerturbation< PairInteraction, ExternalInteraction >::derivative(), McMd::Crosslinker::sample(), and McMd::McPairEnergyAverage::sample().

int McMd::CellList::gridDimension ( int  i) const
inline

Number of cells along axis i.

Parameters
iindex for axis (direction).

Definition at line 466 of file mcMd/neighbor/CellList.h.

int McMd::CellList::totCells ( ) const
inline
int McMd::CellList::nAtom ( ) const

Get total number of atoms in this CellList.

Definition at line 222 of file mcMd/neighbor/CellList.cpp.

Referenced by McMd::PairList::isValid().

int McMd::CellList::atomCapacity ( ) const

Get the maximum allowed atom index + 1.

Definition at line 216 of file mcMd/neighbor/CellList.cpp.

Referenced by McMd::Generator::generate().

bool McMd::CellList::isValid ( int  nAtom = -1) const

Return true if valid, or throw Exception.

If nAtom >= 0, check that number of atoms equals nAtom. If nAtom < 0, skip this check.

Parameters
nAtomnumber of atoms in atom and position arrays.
Returns
true if valid, otherwise throw an Exception.

Definition at line 237 of file mcMd/neighbor/CellList.cpp.

References UTIL_THROW.

Referenced by McMd::PairList::isValid(), and McMd::McSystem::isValid().

Friends And Related Function Documentation

friend class ::CellListTest
friend

Grant access to unit test class.

Definition at line 345 of file mcMd/neighbor/CellList.h.

Member Data Documentation

const int McMd::CellList::MaxNeighbor = 27*Cell::MaxAtomCell
static

Maximum possible number of neighboring atoms.

Use as dimension of the array neighbor passed to getNeighbors() and getCellNeighbors().

Definition at line 69 of file mcMd/neighbor/CellList.h.

Referenced by McMd::GcSliplinkMove::move(), McMd::SliplinkMove::move(), McMd::SliplinkerEnd::move(), McMd::SliplinkerAll::move(), and McMd::Sliplinker::move().


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