ls1-MarDyn
ls1-MarDyn molecular dynamics code
Public Member Functions | Friends | List of all members
LinkedCells Class Reference

Linked Cell Data Structure. More...

#include <LinkedCells.h>

Inheritance diagram for LinkedCells:
ParticleContainer MemoryProfilable

Public Member Functions

 LinkedCells (double bBoxMin[3], double bBoxMax[3], double cutoffRadius)
 initialize the Linked Cell datastructure More...
 
 LinkedCells ()
 Default constructor.
 
 ~LinkedCells ()
 Destructor.
 
void readXML (XMLfileUnits &xmlconfig) override
 Read in XML configuration for LinkedCells and all its included objects. More...
 
int getHaloWidthNumCells () override
 
bool rebuild (double bBoxMin[3], double bBoxMax[3]) override
 rebuild the datastructure More...
 
void update () override
 
void update_via_copies ()
 
void update_via_coloring ()
 
void update_via_traversal ()
 
void update_via_sliced_traversal ()
 
bool addParticle (Molecule &particle, bool inBoxCheckedAlready=false, bool checkWhetherDuplicate=false, const bool &rebuildCaches=false) override
 add a single Molecule to the ParticleContainer. More...
 
void addParticles (std::vector< Molecule > &particles, bool checkWhetherDuplicate=false) override
 adds a whole vector of particles More...
 
void traverseCells (CellProcessor &cellProcessor) override
 calculate the forces between the molecules. More...
 
void traverseNonInnermostCells (CellProcessor &cellProcessor) override
 
void traversePartialInnermostCells (CellProcessor &cellProcessor, unsigned int stage, int stageCount) override
 
unsigned long getNumberOfParticles () override
 
void clear () override
 Delete all molecules in container. More...
 
void deleteOuterParticles () override
 delete all Particles which are not within the bounding box More...
 
double get_halo_L (int index) const override
 gets the width of the halo region in dimension index More...
 
double getCutoff () const override
 
void setCutoff (double rc) override
 
void deleteMolecule (ParticleIterator &moleculeIter, const bool &rebuildCaches) override
 
double getEnergy (ParticlePairsHandler *particlePairsHandler, Molecule *m1, CellProcessor &cellProcessor) override
 
int * getBoxWidthInNumCells ()
 
double * getCellLength () override
 
std::variant< ParticleIterator, SingleCellIterator< ParticleCell > > getMoleculeAtPosition (const double pos[3]) override
 Gets a molecule by its position. More...
 
unsigned long int getCellIndexOfMolecule (Molecule *molecule) const
 Get the index in the cell vector to which this Molecule belongs. More...
 
unsigned long int getCellIndexOfPoint (const double point[3]) const
 Get the index in the cell vector to which the point belongs. More...
 
ParticleCellgetCellReference (int idx)
 
virtual void updateInnerMoleculeCaches () override
 Update the caches of the molecules, that lie in inner cells. The caches of boundary and halo cells is not updated. This method is used for a multi-step scheme of overlapping mpi communication. More...
 
virtual void updateBoundaryAndHaloMoleculeCaches () override
 Update the caches of the molecules, that lie in the boundary or halo cells. The caches of boundary and halo cells is updated, the caches of the inner cells are not updated. This method is used for a multi-step scheme of overlapping mpi communication. More...
 
virtual void updateMoleculeCaches () override
 Update the caches of the molecules. More...
 
ParticleIterator iterator (ParticleIterator::Type t) override
 
RegionParticleIterator regionIterator (const double startRegion[3], const double endRegion[3], ParticleIterator::Type type) override
 
size_t getTotalSize () override
 
void printSubInfo (int offset) override
 
std::string getName () override
 
bool requiresForceExchange () const override
 
unsigned long initCubicGrid (std::array< unsigned long, 3 > numMoleculesPerDimension, std::array< double, 3 > simBoxLength, size_t seed_offset) override
 
std::vector< unsigned long > getParticleCellStatistics () override
 
std::string getConfigurationAsString () override
 
- Public Member Functions inherited from ParticleContainer
 ParticleContainer (double bBoxMin[3], double bBoxMax[3])
 The constructor. More...
 
 ParticleContainer ()
 Default constructor.
 
virtual ~ParticleContainer ()
 The destructor.
 
virtual bool addHaloParticle (Molecule &particle, bool inBoxCheckedAlready=false, bool checkWhetherDuplicate=false, const bool &rebuildCaches=false)
 add a single Molecule to the ParticleContainer, ensures that it is added in the halo. More...
 
virtual double getBoundingBoxMin (int dimension) const
 returns one coordinate of the lower corner of the bounding box More...
 
virtual bool isInBoundingBox (double r[3]) const
 checks, whether given coordinates are within the bounding box More...
 
virtual double getBoundingBoxMax (int dimension) const
 returns one coordinate of the higher corner of the bounding box More...
 
virtual double getSkin () const
 
virtual double * getHaloSize ()
 
virtual std::vector< MoleculegetInvalidParticles ()
 
virtual bool isInvalidParticleReturner ()
 

Friends

class LinkedCellsTest
 
class VTKGridWriter
 

Additional Inherited Members

- Protected Attributes inherited from ParticleContainer
double _boundingBoxMin [3]
 coordinates of the left, lower, front corner of the bounding box
 
double _boundingBoxMax [3]
 coordinates of the right, upper, back corner of the bounding box
 

Detailed Description

Linked Cell Data Structure.

Author
Martin Buchholz

Without any specialized data structure, it needs O(N*N) - where N is the number of particles - time to find all neighbouring pairs of particles. The linked cell data structure is a datastructure which allows to find all neighbouring pairs of particles (neighbouring means particles pairs which have less than a certain distance) in O(N) time. The following picture shows a domain with some particles in it. The blue circle shows the neighbouring area of the red particle The problem is that all particles have to be examined to find those within the circle With the linked cell data structure, the domain is divided into cells (using a regular grid). All particles are placed in those cells. For a given cell, neighbouring cells can easily be calculated, so for a given particle, only the particles from neighbouring cells have to be examined. The following picture illustrates this The spacial domain covered by the linked cells is larger than the bounding box of the domain. This halo region surrounding the phasespace is used for (periodic) boundary conditions and has to be at least as wide as the cutoff radius.
In total, there are three different cell types:

Constructor & Destructor Documentation

◆ LinkedCells()

LinkedCells::LinkedCells ( double  bBoxMin[3],
double  bBoxMax[3],
double  cutoffRadius 
)

initialize the Linked Cell datastructure

The constructor sets the following variables:

  • _cutoffRadius
  • _LJCutoffRadius
  • _haloWidthInNumCells[3]
  • _cellsPerDimension[3]
  • _cellLength[3]
  • _haloBoundingBoxMin and _haloBoundingBoxMax

It resized the cell vector and marks the cells as inner/halo
It fills the array innerCellIndices
It fills the array with forward and backward neighbour indices
The corner parameters for the constructor describe the bounding box of the phasespace which belongs directly to this process, so they correspond to a bounding box including inner + boundary cells but excluding halo cells.
But the corners of this class have to include the halo cells.

Parameters
bBoxMinlower corner of the bounding box of the domain belonging to this container
bBoxMaxhigher corner of the bounding box of the domain belonging to this container
cutoffRadiusdistance for which forces have to be calculated
LJCutoffRadiusdistance for which lennard jones forces have to be calculated
cellsInCutoffRadiusdescribes the width of cells relative to the cutoffRadius:
equal (or larger) to the cutoffRadius divided by the length of a cell as for the number of cells in each dimension only natural numbers are allowed, it can happen that it is not possible to set celllength = cutoffRadius / cellsInCutoffRadius. In that case, the celllength is chosen to be the next larger value so that the sum of the cell lengths in one dimension equals the length of the phasespace Example: phasespacelength=100, cellsInCutoffRadius=2, CutoffRadius=3
==> celllength should be: cutoffRadius/cellsInCutoffRadius = 3/2 = 1.5
==> cellsPerDimension = phasespacelength/celllength = 100/1.5 = 66.67 cells
==> cells have to be larger: cellsPerDimension = phasespacelength/celllength = 100/celllength = 66 cells
==> celllength = 100/66 = 1.5152

Member Function Documentation

◆ addParticle()

bool LinkedCells::addParticle ( Molecule particle,
bool  inBoxCheckedAlready = false,
bool  checkWhetherDuplicate = false,
const bool &  rebuildCaches = false 
)
overridevirtual

add a single Molecule to the ParticleContainer.

Note: a copy of the particle is pushed. Destroying the argument is responsibility of the programmer.

Parameters
particlereference to the particle which has to be added
inBoxCheckedAlready- if true, spare check whether molecule is in bounding box
checkWhetherDuplicate- if true, check whether molecule already exists and don't insert it.
rebuildCachesspecifies, whether the caches should be rebuild
Returns
true if successful, false if particle outside domain

Implements ParticleContainer.

◆ addParticles()

void LinkedCells::addParticles ( std::vector< Molecule > &  particles,
bool  checkWhetherDuplicate = false 
)
overridevirtual

adds a whole vector of particles

Parameters
particlesreference to a vector of pointers to particles

Implements ParticleContainer.

◆ clear()

void LinkedCells::clear ( )
overridevirtual

Delete all molecules in container.

Implements ParticleContainer.

◆ deleteMolecule()

void LinkedCells::deleteMolecule ( ParticleIterator moleculeIter,
const bool &  rebuildCaches 
)
overridevirtual

Implements ParticleContainer.

◆ deleteOuterParticles()

void LinkedCells::deleteOuterParticles ( )
overridevirtual

delete all Particles which are not within the bounding box

Implements ParticleContainer.

◆ get_halo_L()

double LinkedCells::get_halo_L ( int  index) const
overridevirtual

gets the width of the halo region in dimension index

Todo:
remove this method, because a halo_L shouldn't be necessary for every ParticleContainer e.g. replace it by the cutoff-radius

Implements ParticleContainer.

◆ getCellIndexOfMolecule()

unsigned long int LinkedCells::getCellIndexOfMolecule ( Molecule molecule) const

Get the index in the cell vector to which this Molecule belongs.

each spatial position within the bounding box of the linked cells belongs unambiguously to one cell.
This method determines for a given Molecule the corresponding cell and returns the index of that cell in the cell vector.
If the molecule is not inside the bounding box, an error is printed

◆ getCellIndexOfPoint()

unsigned long int LinkedCells::getCellIndexOfPoint ( const double  point[3]) const

Get the index in the cell vector to which the point belongs.

each spatial position within the bounding box of the linked cells belongs unambiguously to one cell.
This method determines for a given point the corresponding cell and returns the index of that cell in the cell vector.
If the point is not inside the bounding box, an error is printed

◆ getCellLength()

double* LinkedCells::getCellLength ( )
inlineoverridevirtual

Implements ParticleContainer.

◆ getConfigurationAsString()

string LinkedCells::getConfigurationAsString ( )
overridevirtual

Return a string representation of the algorithmic configuration of the container. Only used for logging / output.

Note
: Formatting rules:
  • The whole configuration should be enclosed in curly brackets.
  • Different elements should be separated by a comma sourrounded by spaces: " , ".
  • Every element should be in the form "key: value"
Returns

Implements ParticleContainer.

◆ getCutoff()

double LinkedCells::getCutoff ( ) const
inlineoverridevirtual

Implements ParticleContainer.

◆ getEnergy()

double LinkedCells::getEnergy ( ParticlePairsHandler particlePairsHandler,
Molecule m1,
CellProcessor cellProcessor 
)
overridevirtual

Implements ParticleContainer.

◆ getHaloWidthNumCells()

int LinkedCells::getHaloWidthNumCells ( )
inlineoverridevirtual

Reimplemented from ParticleContainer.

◆ getMoleculeAtPosition()

std::variant< ParticleIterator, SingleCellIterator< ParticleCell > > LinkedCells::getMoleculeAtPosition ( const double  pos[3])
overridevirtual

Gets a molecule by its position.

Parameters
posMolecule position
resultMolecule will be returned by this pointer if found
Returns
Molecule was found?

Implements ParticleContainer.

◆ getName()

std::string LinkedCells::getName ( )
overridevirtual

Implements MemoryProfilable.

◆ getNumberOfParticles()

unsigned long LinkedCells::getNumberOfParticles ( )
overridevirtual
Returns
the number of particles stored in the Linked Cells

Implements ParticleContainer.

◆ getParticleCellStatistics()

std::vector< unsigned long > LinkedCells::getParticleCellStatistics ( )
overridevirtual

Get a statistics of the cells found in this container.

Returns
A vector, where particleCellStatistics[partCount] represents the number of cells with partCount particles.

Reimplemented from ParticleContainer.

◆ getTotalSize()

size_t LinkedCells::getTotalSize ( )
overridevirtual

Implements MemoryProfilable.

◆ initCubicGrid()

unsigned long LinkedCells::initCubicGrid ( std::array< unsigned long, 3 >  numMoleculesPerDimension,
std::array< double, 3 >  simBoxLength,
size_t  seed_offset 
)
overridevirtual

Implements ParticleContainer.

◆ iterator()

ParticleIterator LinkedCells::iterator ( ParticleIterator::Type  t)
inlineoverridevirtual

Implements ParticleContainer.

◆ printSubInfo()

void LinkedCells::printSubInfo ( int  offset)
overridevirtual

Implements MemoryProfilable.

◆ readXML()

void LinkedCells::readXML ( XMLfileUnits xmlconfig)
overridevirtual

Read in XML configuration for LinkedCells and all its included objects.

The following xml object structure is handled by this method:

<datastructure type="LinkedCells">
<cellsInCutoffRadius>INTEGER</cellsInCutoffRadius>
<!-- from TraversalTuner: -->
<!-- select traversal algorithm
possible values are:
- original
- c04
- c08 (default for >1 threads)
- c08es (eight-shell)
- quicksched
- sliced (default for <2 threads)
- hs (half shell method)
- mp (mid point method)
- nt (neutral territory method)
-->
<traversalSelector>c08</traversalSelector>
<!-- override default block size (2x2x2) for quicksched tasks -->
<traversalData type="quicksched">
<taskBlockSize>
<lx>2</lx>
<ly>2</ly>
<lz>2</lz>
</taskBlockSize>
</traversalData>
</datastructure>

Implements ParticleContainer.

◆ rebuild()

bool LinkedCells::rebuild ( double  bBoxMin[3],
double  bBoxMax[3] 
)
overridevirtual

rebuild the datastructure

Load-balancing decompositions change the position and size of the local region during runtime. Therefore, the datastructure needs to be rebuild completely. This method basically does what the constructor does as well, with the difference, that there are already particles stored, and particles which don't belong to the new region have to be deleted after rebuild @parameter bBoxMin minimum of the box @parameter bBoxMax maximum of the box

Reimplemented from ParticleContainer.

◆ regionIterator()

RegionParticleIterator LinkedCells::regionIterator ( const double  startRegion[3],
const double  endRegion[3],
ParticleIterator::Type  type 
)
overridevirtual

Implements ParticleContainer.

◆ requiresForceExchange()

bool LinkedCells::requiresForceExchange ( ) const
overridevirtual

Reimplemented from ParticleContainer.

◆ setCutoff()

void LinkedCells::setCutoff ( double  rc)
inlineoverridevirtual

set the cutoff

Parameters
rc

Reimplemented from ParticleContainer.

◆ traverseCells()

void LinkedCells::traverseCells ( CellProcessor cellProcessor)
overridevirtual

calculate the forces between the molecules.

Only molecules with a distance not larger than the cutoff radius are to be used.
Only forces on the Molecules which are in the inner and boundary region have to be calculated Newton's third law should be used for faster computation:

  • a loop over all inner cells calculates the forces with all forward neighbour cells all forward cells have to be used, as none of them can be halo or outside
  • a loop over the boundary cells first calculates forces with all forward cells and all backward cells. Here it has to be checked whether the neighbour cell is halo or not. If it is Halo, the force is calculated, if it isn't, the force is not calculated, because the same pair of cells has already been processed in one of the other loops.
    Parameters
    particlePairsHandlerspecified concrete action to be done for each pair

Implements ParticleContainer.

◆ traverseNonInnermostCells()

void LinkedCells::traverseNonInnermostCells ( CellProcessor cellProcessor)
overridevirtual

Implements ParticleContainer.

◆ traversePartialInnermostCells()

void LinkedCells::traversePartialInnermostCells ( CellProcessor cellProcessor,
unsigned int  stage,
int  stageCount 
)
overridevirtual

Implements ParticleContainer.

◆ update()

void LinkedCells::update ( )
overridevirtual

Pointers to the particles are put into cells depending on the spacial position of the particles. Before the call of this method, this distribution might have become invalid. To ensure, that all Particles (pointers to them) are put into the corresponding cells, first all cells are cleared and then filled again depending on the spacial position of the molecules. After the update, exactly one pointer for each particle in this ParticleContainer is it's corresponding cell.

Implements ParticleContainer.

◆ updateBoundaryAndHaloMoleculeCaches()

void LinkedCells::updateBoundaryAndHaloMoleculeCaches ( )
overridevirtual

Update the caches of the molecules, that lie in the boundary or halo cells. The caches of boundary and halo cells is updated, the caches of the inner cells are not updated. This method is used for a multi-step scheme of overlapping mpi communication.

Implements ParticleContainer.

◆ updateInnerMoleculeCaches()

void LinkedCells::updateInnerMoleculeCaches ( )
overridevirtual

Update the caches of the molecules, that lie in inner cells. The caches of boundary and halo cells is not updated. This method is used for a multi-step scheme of overlapping mpi communication.

Implements ParticleContainer.

◆ updateMoleculeCaches()

void LinkedCells::updateMoleculeCaches ( )
overridevirtual

Update the caches of the molecules.

Implements ParticleContainer.


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