ls1-MarDyn
ls1-MarDyn molecular dynamics code
Public Member Functions | Protected Attributes | List of all members
ParticleContainer Class Referenceabstract

This Interface is used to get access to particles and pairs of particles. More...

#include <ParticleContainer.h>

Inheritance diagram for ParticleContainer:
MemoryProfilable AutoPasContainer LinkedCells ParticleContainerToBasisWrapper

Public Member Functions

 ParticleContainer (double bBoxMin[3], double bBoxMax[3])
 The constructor. More...
 
 ParticleContainer ()
 Default constructor.
 
virtual ~ParticleContainer ()
 The destructor.
 
virtual void readXML (XMLfileUnits &xmlconfig)=0
 
virtual bool rebuild (double bBoxMin[3], double bBoxMax[3])
 rebuild the datastructure More...
 
virtual void update ()=0
 do necessary updates resulting from changed particle positions More...
 
virtual bool addParticle (Molecule &particle, bool inBoxCheckedAlready=false, bool checkWhetherDuplicate=false, const bool &rebuildCaches=false)=0
 add a single Molecule to the ParticleContainer. More...
 
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 void addParticles (std::vector< Molecule > &particles, bool checkWhetherDuplicate=false)=0
 adds a whole vector of particles More...
 
virtual void traverseCells (CellProcessor &cellProcessor)=0
 traverse pairs which are close to each other More...
 
virtual void traverseNonInnermostCells (CellProcessor &cellProcessor)=0
 
virtual void traversePartialInnermostCells (CellProcessor &cellProcessor, unsigned int stage, int stageCount)=0
 
virtual ParticleIterator iterator (ParticleIterator::Type t)=0
 
virtual RegionParticleIterator regionIterator (const double startCorner[3], const double endCorner[3], ParticleIterator::Type t)=0
 
virtual unsigned long getNumberOfParticles ()=0
 
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 int getHaloWidthNumCells ()
 
virtual double getBoundingBoxMax (int dimension) const
 returns one coordinate of the higher corner of the bounding box More...
 
virtual void clear ()=0
 Delete all molecules in container. More...
 
virtual void deleteOuterParticles ()=0
 delete all Particles which are not within the bounding box More...
 
virtual double get_halo_L (int index) const =0
 returns the width of the halo stripe (for the given dimension index) More...
 
virtual double getCutoff () const =0
 
virtual double getSkin () const
 
virtual void deleteMolecule (ParticleIterator &moleculeIter, const bool &rebuildCaches)=0
 
virtual double getEnergy (ParticlePairsHandler *particlePairsHandler, Molecule *m1, CellProcessor &cellProcessor)=0
 
virtual void updateInnerMoleculeCaches ()=0
 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 ()=0
 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 ()=0
 Update the caches of the molecules. More...
 
virtual std::variant< ParticleIterator, SingleCellIterator< ParticleCell > > getMoleculeAtPosition (const double pos[3])=0
 Gets a molecule by its position. More...
 
virtual bool requiresForceExchange () const
 
virtual unsigned long initCubicGrid (std::array< unsigned long, 3 > numMoleculesPerDimension, std::array< double, 3 > simBoxLength, size_t seed_offset)=0
 
virtual double * getCellLength ()=0
 
virtual double * getHaloSize ()
 
virtual std::vector< unsigned long > getParticleCellStatistics ()
 
virtual void setCutoff (double rc)
 
virtual std::vector< MoleculegetInvalidParticles ()
 
virtual bool isInvalidParticleReturner ()
 
virtual std::string getConfigurationAsString ()=0
 
- Public Member Functions inherited from MemoryProfilable
virtual size_t getTotalSize ()=0
 
virtual void printSubInfo (int offset)=0
 
virtual std::string getName ()=0
 

Protected Attributes

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

This Interface is used to get access to particles and pairs of particles.

Author
Martin Buchholz

A particleContainer is used to store Particles with a short-range potential in a way that the access to pairs of neighbouring particles is efficient. Neighbouring particles are particles which have a distance shorter than a given cutoff radius. A common task when using a PariticleContainer is to do something for all particles. This can be done using the methods begin, next and end, e.g.:

ParticleContainer* partContPtr;
Molecule* partPtr;
for(partPtr = partContPtr->begin(); partPtr != partContPtr->end(); partPtr = partContPtr->next()){
partPtr->doSomething();
}
FullMolecule modeled as LJ sphere with point polarities.
Definition: FullMolecule.h:18
This Interface is used to get access to particles and pairs of particles.
Definition: ParticleContainer.h:69

Particles stored in this container can either belong to the local process or they can be duplicates (from neighbouring processes or periodic boundary). As the simulated regions are often cuboids, a bounding box is defined through two opposing corners (_boundingBoxMin[3] and _boundingBoxMax[3]). Particles inside this bounding box belong to this process, those outside don't. An exception to this is when particles are moved in a time step. It has to be ensured that particles which leave the bounding box are properly handled.

For non-cuboid regions, the bounding box still has to be defined as it gives an approximation for the region that is covered by the ParticleContainer.

This interface doesn't implement the datastructure, it just tells which methods a class implementing this kind of datastructure has to provide to be used by the framework. Such a class should be implemented as a subclass of this class.

Constructor & Destructor Documentation

◆ ParticleContainer()

ParticleContainer::ParticleContainer ( double  bBoxMin[3],
double  bBoxMax[3] 
)

The constructor.

Parameters
bBoxMincoordinates of the lowest (in all coordinates) corner of the bounding box
bBoxMaxcoordinates of the highest (in all coordinates) corner of the bounding box

Member Function Documentation

◆ addHaloParticle()

bool ParticleContainer::addHaloParticle ( Molecule particle,
bool  inBoxCheckedAlready = false,
bool  checkWhetherDuplicate = false,
const bool &  rebuildCaches = false 
)
virtual

add a single Molecule to the ParticleContainer, ensures that it is added in the halo.

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

Reimplemented in AutoPasContainer.

◆ addParticle()

virtual bool ParticleContainer::addParticle ( Molecule particle,
bool  inBoxCheckedAlready = false,
bool  checkWhetherDuplicate = false,
const bool &  rebuildCaches = false 
)
pure virtual

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

Implemented in AutoPasContainer, LinkedCells, and ParticleContainerToBasisWrapper.

◆ addParticles()

virtual void ParticleContainer::addParticles ( std::vector< Molecule > &  particles,
bool  checkWhetherDuplicate = false 
)
pure virtual

adds a whole vector of particles

Parameters
particlesreference to a vector of pointers to particles

Implemented in AutoPasContainer, LinkedCells, and ParticleContainerToBasisWrapper.

◆ clear()

virtual void ParticleContainer::clear ( )
pure virtual

Delete all molecules in container.

Implemented in AutoPasContainer, LinkedCells, and ParticleContainerToBasisWrapper.

◆ deleteOuterParticles()

virtual void ParticleContainer::deleteOuterParticles ( )
pure virtual

delete all Particles which are not within the bounding box

Implemented in AutoPasContainer, LinkedCells, and ParticleContainerToBasisWrapper.

◆ get_halo_L()

virtual double ParticleContainer::get_halo_L ( int  index) const
pure virtual

returns the width of the halo stripe (for the given 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

Implemented in AutoPasContainer, LinkedCells, and ParticleContainerToBasisWrapper.

◆ getBoundingBoxMax()

double ParticleContainer::getBoundingBoxMax ( int  dimension) const
virtual

returns one coordinate of the higher corner of the bounding box

Parameters
dimensionthe coordinate which should be returned

Reimplemented in ParticleContainerToBasisWrapper.

◆ getBoundingBoxMin()

double ParticleContainer::getBoundingBoxMin ( int  dimension) const
virtual

returns one coordinate of the lower corner of the bounding box

Parameters
dimensionthe coordinate which should be returned

Reimplemented in ParticleContainerToBasisWrapper.

◆ getConfigurationAsString()

virtual std::string ParticleContainer::getConfigurationAsString ( )
pure virtual

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

Implemented in AutoPasContainer, LinkedCells, and ParticleContainerToBasisWrapper.

◆ getHaloSize()

virtual double* ParticleContainer::getHaloSize ( )
inlinevirtual

Return the size of the halo

Returns

Reimplemented in AutoPasContainer.

◆ getMoleculeAtPosition()

virtual std::variant<ParticleIterator, SingleCellIterator<ParticleCell> > ParticleContainer::getMoleculeAtPosition ( const double  pos[3])
pure virtual

Gets a molecule by its position.

Parameters
posMolecule position
Returns
Iterator to the molecule. The iterator is invalid if no molecule was found.

Implemented in AutoPasContainer, LinkedCells, and ParticleContainerToBasisWrapper.

◆ getNumberOfParticles()

virtual unsigned long ParticleContainer::getNumberOfParticles ( )
pure virtual
Returns
the number of particles stored in this container

This number may include particles which are outside of the bounding box

Implemented in AutoPasContainer, LinkedCells, and ParticleContainerToBasisWrapper.

◆ getParticleCellStatistics()

virtual std::vector<unsigned long> ParticleContainer::getParticleCellStatistics ( )
inlinevirtual

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 in LinkedCells.

◆ isInBoundingBox()

bool ParticleContainer::isInBoundingBox ( double  r[3]) const
virtual

checks, whether given coordinates are within the bounding box

Parameters
rthe coordinates to check
Returns
true if coordinates within bounding box, false otherwise

Reimplemented in ParticleContainerToBasisWrapper.

◆ readXML()

virtual void ParticleContainer::readXML ( XMLfileUnits xmlconfig)
pure virtual

Implemented in AutoPasContainer, and LinkedCells.

◆ rebuild()

bool ParticleContainer::rebuild ( double  bBoxMin[3],
double  bBoxMax[3] 
)
virtual

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 in AutoPasContainer, and LinkedCells.

◆ setCutoff()

virtual void ParticleContainer::setCutoff ( double  rc)
inlinevirtual

set the cutoff

Parameters
rc

Reimplemented in AutoPasContainer, and LinkedCells.

◆ traverseCells()

virtual void ParticleContainer::traverseCells ( CellProcessor cellProcessor)
pure virtual

traverse pairs which are close to each other

Only interactions between particles which have a distance which is not larger than the cutoff radius are to be considered.
This method has to be implemented in derived classes Precondition: All Particles of the process + halo molecules are stored Task: Run over all pairs (Each pair exactely once!) of particles (within cutoffradius) Important: Some pairs might be "duplicated": All pairs which cross the boundary occur twice (second time at the periodic image). Those pairs are from the point of view of the datastructure two different pairs, but they both times connect the same particles. For a pair which occurs twice, it has to be made sure, that one gets the status "original pair" and the other one "duplicated pair". For each pair found, there is an action executed, but it is a different action for original and duplicated pairs. Details about how to handle pairs can be found in the documentation for the class ParticlePairsHandler

Parameters
particlePairsHandlerspecified concrete action to be done for each pair

Implemented in AutoPasContainer, LinkedCells, and ParticleContainerToBasisWrapper.

◆ update()

virtual void ParticleContainer::update ( )
pure virtual

do necessary updates resulting from changed particle positions

For some implementations of the interface ParticleContainer, the place where particles are stored might e.g. depend on the spacial position of the particle. So when some externel method (e.g. Leap-Frog) changes the spacial position of a particle, the representation within the particleContainer becomes invalid. This method restores a valid representation.

Implemented in AutoPasContainer, LinkedCells, and ParticleContainerToBasisWrapper.

◆ updateBoundaryAndHaloMoleculeCaches()

virtual void ParticleContainer::updateBoundaryAndHaloMoleculeCaches ( )
pure virtual

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.

Implemented in AutoPasContainer, LinkedCells, and ParticleContainerToBasisWrapper.

◆ updateInnerMoleculeCaches()

virtual void ParticleContainer::updateInnerMoleculeCaches ( )
pure virtual

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.

Implemented in AutoPasContainer, LinkedCells, and ParticleContainerToBasisWrapper.

◆ updateMoleculeCaches()

virtual void ParticleContainer::updateMoleculeCaches ( )
pure virtual

Update the caches of the molecules.

Implemented in AutoPasContainer, LinkedCells, and ParticleContainerToBasisWrapper.


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