Molecular Dynamics Simulation  1.0
Public Types | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes | Friends | List of all members
LinkedCellsContainer Class Reference

Extension of the ParticleContainer class using a linked cells data structure for improved performance over the direct sum approach. More...

#include <LinkedCellsContainer.h>

Inheritance diagram for LinkedCellsContainer:
Inheritance graph
Collaboration diagram for LinkedCellsContainer:
Collaboration graph

Public Types

enum class  BoundaryCondition { OUTFLOW , REFLECTIVE , PERIODIC }
 Boundary type enum for labeling the sides of the domain. More...
 
enum class  BoundarySide {
  LEFT , RIGHT , BOTTOM , TOP ,
  BACK , FRONT
}
 Boundary side enum for labeling the sides of the domain. More...
 

Public Member Functions

 LinkedCellsContainer (const std::array< double, 3 > &domain_size, double cutoff_radius, const std::array< BoundaryCondition, 6 > &boundary_types={BoundaryCondition::OUTFLOW, BoundaryCondition::OUTFLOW, BoundaryCondition::OUTFLOW, BoundaryCondition::OUTFLOW, BoundaryCondition::OUTFLOW, BoundaryCondition::OUTFLOW}, int n=0)
 Construct a new Linked Cells Particle Container object. More...
 
void addParticle (const Particle &p) override
 Adds a particle to the container. More...
 
void addParticle (Particle &&p) override
 Adds a particle to the container. More...
 
void prepareForceCalculation () override
 Prepares everything for the force calculations (must be called before applySimpleForces and applyPairwiseForces) More...
 
void applySimpleForces (const std::vector< std::shared_ptr< SimpleForceSource >> &simple_force_sources) override
 Applies the given simple force sources to the particles. More...
 
void applyPairwiseForces (const std::vector< std::shared_ptr< PairwiseForceSource >> &force_sources) override
 Applies the given force sources to the particles. More...
 
void applyTargettedForces (const std::vector< std::shared_ptr< TargettedForceSource >> &targetted_force_sources, double curr_simulation_time) override
 Applies the given targetted force sources to the particles. More...
 
void reserve (size_t n) override
 Reserves space for n particles. This is useful if the number of particles is known in advance and prevents reallocation of memory for the internal dynamic array of particles, when inserting new particles. More...
 
size_t size () const override
 Returns the number of particles in the container. More...
 
Particleoperator[] (int i) override
 Overload of the [] operator to access the particles in the container. More...
 
std::vector< Particle >::iterator begin () override
 The begin iterator for the internal data structure. More...
 
std::vector< Particle >::iterator end () override
 The end iterator for the internal data structure. More...
 
std::vector< Particle >::const_iterator begin () const override
 The begin const iterator for the internal data structure. More...
 
std::vector< Particle >::const_iterator end () const override
 The end const iterator for the internal data structure. More...
 
const std::vector< Particle > & getParticles () const override
 Returns a vector of all particles in the container. More...
 
const std::array< double, 3 > & getDomainSize () const
 Returns the domain size. More...
 
double getCutoffRadius () const
 Returns the cutoff radius. More...
 
const std::vector< Cell > & getCells ()
 Returns the cells. More...
 
const std::vector< Cell * > & getBoundaryCells () const
 Returns the pointers of the boundary cells. More...
 
const std::array< double, 3 > & getCellSize () const
 Returns the cell size. More...
 
const std::array< int, 3 > & getDomainNumCells () const
 Returns the number of cells in each dimension. More...
 
int cellCoordToCellIndex (int cx, int cy, int cz) const
 Maps the cell coordinates to the corresponding index in the internal cell vector. More...
 
CellparticlePosToCell (const std::array< double, 3 > &pos)
 Maps the particle position to the corresponding cell index in the internal cell vector. More...
 
CellparticlePosToCell (double x, double y, double z)
 Maps the particle position to the corresponding cell index in the internal cell vector. More...
 
- Public Member Functions inherited from ParticleContainer
virtual ~ParticleContainer ()=default
 Virtual destructor for correct deconstruction of inheriting classes. More...
 

Static Public Member Functions

static std::string boundaryConditionToString (const BoundaryCondition &bc)
 Returns a string description of a boundary condition. More...
 

Private Member Functions

void initCells ()
 Populates the cell vector and sets the cells types. More...
 
void initCellNeighbourReferences ()
 Sets the neighbour references for each cell in the cell vector. More...
 
void initIterationOrders ()
 Sets the iteration orders for parallel execution. More...
 
void updateCellsParticleReferences ()
 Updates the particle references in the cells. This is necessary after a reallocation of the internal particle vector. More...
 
void deleteHaloParticles ()
 Removes all particles in the halo cells from the particles vector. ATTENTION: A particle reference update must be triggered manually after this operation! More...
 

Private Attributes

std::vector< std::vector< Cell * > > iteration_orders
 
std::vector< Particleparticles
 Internal data structure for the particles. More...
 
std::array< double, 3 > domain_size
 Domain size in each dimension. More...
 
double cutoff_radius
 Cutoff radius for the force calculation. More...
 
std::array< BoundaryCondition, 6 > boundary_types
 The boundary types for each side of the domain (order in array: left, right, bottom, top, back, front) More...
 
std::array< double, 3 > cell_size
 Cell size in each dimension. More...
 
std::array< int, 3 > domain_num_cells
 Number of cells in each dimension. More...
 
std::vector< Cellcells
 Internal data structure for the cells. More...
 
std::vector< Cell * > domain_cell_references
 References to the domain cells. More...
 
std::vector< Cell * > boundary_cell_references
 References to the boundary cells. More...
 
std::vector< Cell * > halo_cell_references
 References to the halo cells. More...
 
std::vector< Cell * > left_boundary_cell_references
 References to the boundary cells on the left (x = 0) More...
 
std::vector< Cell * > right_boundary_cell_references
 References to the boundary cells on the right (x = domain_num_cells[0]-1) More...
 
std::vector< Cell * > bottom_boundary_cell_references
 References to the boundary cells on the bottom (y = 0) More...
 
std::vector< Cell * > top_boundary_cell_references
 References to the boundary cells on the top (y = domain_num_cells[1]-1) More...
 
std::vector< Cell * > back_boundary_cell_references
 References to the boundary cells on the back (z = 0) More...
 
std::vector< Cell * > front_boundary_cell_references
 References to the boundary cells on the front (z = domain_num_cells[2]-1) More...
 
std::vector< Cell * > left_halo_cell_references
 References to the halo cells on the left (x = -1) More...
 
std::vector< Cell * > right_halo_cell_references
 References to the halo cells on the right (x = domain_num_cells[0]) More...
 
std::vector< Cell * > bottom_halo_cell_references
 References to the halo cells on the bottom (y = -1) More...
 
std::vector< Cell * > top_halo_cell_references
 References to the halo cells on the top (y = domain_num_cells[1]) More...
 
std::vector< Cell * > back_halo_cell_references
 References to the halo cells on the back (z = -1) More...
 
std::vector< Cell * > front_halo_cell_references
 References to the halo cells on the front (z = domain_num_cells[2]) More...
 

Friends

class ReflectiveBoundaryType
 Friend class to allow access to the internal data structures. More...
 
class OutflowBoundaryType
 Friend class to allow access to the internal data structures. More...
 
class PeriodicBoundaryType
 Friend class to allow access to the internal data structures. More...
 

Detailed Description

Extension of the ParticleContainer class using a linked cells data structure for improved performance over the direct sum approach.

Definition at line 14 of file LinkedCellsContainer.h.

Member Enumeration Documentation

◆ BoundaryCondition

Boundary type enum for labeling the sides of the domain.

Enumerator
OUTFLOW 
REFLECTIVE 
PERIODIC 

Definition at line 21 of file LinkedCellsContainer.h.

21 { OUTFLOW, REFLECTIVE, PERIODIC };

◆ BoundarySide

Boundary side enum for labeling the sides of the domain.

Enumerator
LEFT 
RIGHT 
BOTTOM 
TOP 
BACK 
FRONT 

Definition at line 26 of file LinkedCellsContainer.h.

26 { LEFT, RIGHT, BOTTOM, TOP, BACK, FRONT };

Constructor & Destructor Documentation

◆ LinkedCellsContainer()

LinkedCellsContainer::LinkedCellsContainer ( const std::array< double, 3 > &  domain_size,
double  cutoff_radius,
const std::array< BoundaryCondition, 6 > &  boundary_types = {BoundaryCondition::OUTFLOWBoundaryCondition::OUTFLOWBoundaryCondition::OUTFLOWBoundaryCondition::OUTFLOWBoundaryCondition::OUTFLOWBoundaryCondition::OUTFLOW},
int  n = 0 
)

Construct a new Linked Cells Particle Container object.

Parameters
domain_sizethe size of the domain
cutoff_radiusthe cutoff radius for the force calculation
boundary_typesthe boundary types for each side of the domain (order in array: left, right, bottom, top, back, front)
nthe expected number of particles (for preallocation of memory)

Constructs a new Linked Cells Particle Container object using the specified domain size and cutoff radius. The expected number of particles is used to preallocate memory for the particle vector, which is highly recommended for known amounts of particles (improves performance as no std::vector resize is needed on inserts). The origin (left lower front corner of the boundary) is assumed to be at (0, 0, 0). The cutoff radius is used to determine the number of cells in each dimension, where the cell size is bigger or equal then the cutoff radius but adjusted to divide the domain size evenly into a whole number of cells. Ideally the domain size is a multiple of the cutoff and cube shaped. Internally, the domain is extended by one cell in each direction to allow for so called 'ghost' or 'halo' cells that are used for boundary condition handling. Therefore the valid cell coordinates range from -1 to domain_num_cells[i] in each dimension (i = 0 -> x; i = 1 -> y; i = 2 -> z).

Definition at line 18 of file LinkedCellsContainer.cpp.

20  : domain_size(_domain_size), cutoff_radius(_cutoff_radius), boundary_types(_boundary_types) {
21  domain_num_cells = {std::max(static_cast<int>(std::floor(_domain_size[0] / cutoff_radius)), 1),
22  std::max(static_cast<int>(std::floor(_domain_size[1] / cutoff_radius)), 1),
23  std::max(static_cast<int>(std::floor(_domain_size[2] / cutoff_radius)), 1)};
24 
25  cell_size = {_domain_size[0] / domain_num_cells[0], _domain_size[1] / domain_num_cells[1], _domain_size[2] / domain_num_cells[2]};
26 
27  // reserve the memory for the cells
28  cells.reserve((domain_num_cells[0] + 2) * (domain_num_cells[1] + 2) * (domain_num_cells[2] + 2));
29 
30  // create the cells with the correct cell-type and add them to the cells vector and the corresponding cell reference vector
31  initCells();
32 
33  // add the neighbour references to the cells
35 
36  // create the iteration orders
38 
39  // reserve the memory for the particles to prevent reallocation during insertion
40  particles.reserve(_n);
41 
42  Logger::logger->debug("Created LinkedCellsContainer with boundaries [{}, {}, {}] and cutoff radius {}", domain_size[0], domain_size[1],
44  Logger::logger->debug("Created LinkedCellsContainer with {} domain cells (of which {} are at the boundary) and {} halo cells",
46  Logger::logger->debug("Cells per dimension: [{}, {}, {}]", domain_num_cells[0], domain_num_cells[1], domain_num_cells[2]);
47  Logger::logger->debug("Calculated cell size: [{}, {}, {}]", cell_size[0], cell_size[1], cell_size[2]);
48 }
std::vector< Cell * > domain_cell_references
References to the domain cells.
std::vector< Cell * > boundary_cell_references
References to the boundary cells.
void initIterationOrders()
Sets the iteration orders for parallel execution.
std::array< double, 3 > cell_size
Cell size in each dimension.
std::array< double, 3 > domain_size
Domain size in each dimension.
std::array< int, 3 > domain_num_cells
Number of cells in each dimension.
std::vector< Cell > cells
Internal data structure for the cells.
void initCells()
Populates the cell vector and sets the cells types.
void initCellNeighbourReferences()
Sets the neighbour references for each cell in the cell vector.
std::vector< Cell * > halo_cell_references
References to the halo cells.
double cutoff_radius
Cutoff radius for the force calculation.
std::array< BoundaryCondition, 6 > boundary_types
The boundary types for each side of the domain (order in array: left, right, bottom,...
std::vector< Particle > particles
Internal data structure for the particles.
static std::shared_ptr< spdlog::logger > logger
Publically accessible shared pointer to the logger.
Definition: Logger.h:35

Member Function Documentation

◆ addParticle() [1/2]

void LinkedCellsContainer::addParticle ( const Particle p)
overridevirtual

Adds a particle to the container.

Parameters
pParticle to be added

Adds a particle to the container and correctly inserts it into the cell structure.

Implements ParticleContainer.

Definition at line 50 of file LinkedCellsContainer.cpp.

50  {
51  Cell* cell = particlePosToCell(p.getX());
52 
53  if (cell == nullptr) {
54  Logger::logger->error("Particle to insert is out of bounds, position: [{}, {}, {}]", p.getX()[0], p.getX()[1], p.getX()[2]);
55  throw std::runtime_error("Attempted to insert particle out of bounds");
56  }
57  // if (cell->getCellType() == Cell::CellType::HALO) {
58  // Logger::logger->warn("Particle to insert is in halo cell. Position: [{}, {}, {}]", p.getX()[0], p.getX()[1], p.getX()[2]);
59  // throw std::runtime_error("Attempted to insert particle into halo cell");
60  // }
61 
62  size_t old_capacity = particles.capacity();
63  particles.push_back(p);
64 
65  if (old_capacity != particles.capacity()) {
67  } else {
68  cell->addParticleReference(&particles.back());
69  }
70 }
Class representing a cell in the linked cells algorithm.
Definition: Cell.h:12
void addParticleReference(Particle *p)
Adds a particle reference to the cell.
Definition: Cell.cpp:11
void updateCellsParticleReferences()
Updates the particle references in the cells. This is necessary after a reallocation of the internal ...
Cell * particlePosToCell(const std::array< double, 3 > &pos)
Maps the particle position to the corresponding cell index in the internal cell vector.
const std::array< double, 3 > & getX() const
Gets the position of the particle.
Definition: Particle.h:157

◆ addParticle() [2/2]

void LinkedCellsContainer::addParticle ( Particle &&  p)
overridevirtual

Adds a particle to the container.

Parameters
pParticle to be added

Adds a particle to the container and correctly inserts it into the cell structure.

Implements ParticleContainer.

Definition at line 72 of file LinkedCellsContainer.cpp.

72  {
73  Cell* cell = particlePosToCell(p.getX());
74 
75  if (cell == nullptr) {
76  Logger::logger->error("Particle to insert is outside of cells. Position: [{}, {}, {}]", p.getX()[0], p.getX()[1], p.getX()[2]);
77  throw std::runtime_error("Attempted to insert particle out of bounds");
78  }
79  // if (cell->getCellType() == Cell::CellType::HALO) {
80  // Logger::logger->warn("Particle to insert is in halo cell. Position: [{}, {}, {}]", p.getX()[0], p.getX()[1], p.getX()[2]);
81  // throw std::runtime_error("Attempted to insert particle into halo cell");
82  // }
83 
84  size_t old_capacity = particles.capacity();
85  particles.push_back(std::move(p));
86 
87  if (old_capacity != particles.capacity()) {
89  } else {
90  cell->addParticleReference(&particles.back());
91  }
92 }

◆ applyPairwiseForces()

void LinkedCellsContainer::applyPairwiseForces ( const std::vector< std::shared_ptr< PairwiseForceSource >> &  force_sources)
overridevirtual

Applies the given force sources to the particles.

Parameters
force_sourcesList of force sources to be applied

Applies the given force sources to the particles in the container. Uses newton's third law to calculate the forces between the particles in an optimized way. Additionally to the functionality of the ParticleContainer class, this method uses the internal cell structure to reduce the number of force calculations necessary, depending on the cutoff radius.

Implements ParticleContainer.

Definition at line 117 of file LinkedCellsContainer.cpp.

117  {
118  // apply the boundary conditions
120  size_t original_particle_length = PeriodicBoundaryType::applyBoundaryConditions(*this);
121 
122  for (auto& it_order : iteration_orders) {
123 #ifdef _OPENMP
124 #pragma omp parallel for schedule(dynamic)
125 #endif
126  for (Cell* cell : it_order) {
127  if (cell->getParticleReferences().empty()) continue;
128 
129  // Calculate cell internal forces
130  for (auto it1 = cell->getParticleReferences().begin(); it1 != cell->getParticleReferences().end(); ++it1) {
131  Particle* p = *it1;
132  // calculate the forces between the particle and the particles in the same cell
133  // uses direct sum with newtons third law
134  for (auto it2 = (it1 + 1); it2 != cell->getParticleReferences().end(); ++it2) {
135  if (p->isLocked() && (*it2)->isLocked()) continue;
136  Particle* q = *it2;
137  std::array<double, 3> total_force{0, 0, 0};
138  for (auto& force : force_sources) {
139  total_force = total_force + force->calculateForce(*p, *q);
140  }
141  p->addF(total_force);
142  q->subF(total_force);
143  }
144  }
145 
146  // Calculate cell boundary forces
147  for (Cell* neighbour_cell : cell->getNeighbourReferences()) {
148  if (cell < neighbour_cell) continue;
149  if (neighbour_cell->getParticleReferences().empty()) continue;
150 
151  for (Particle* p : cell->getParticleReferences()) {
152  for (Particle* neighbour_particle : neighbour_cell->getParticleReferences()) {
153  if (p->isLocked() && neighbour_particle->isLocked()) continue;
154  if (ArrayUtils::L2NormSquared(p->getX() - neighbour_particle->getX()) > cutoff_radius * cutoff_radius) continue;
155 
156  auto total_force = std::array<double, 3>{0, 0, 0};
157  for (const auto& force_source : force_sources) {
158  total_force = total_force + force_source->calculateForce(*p, *neighbour_particle);
159  }
160  p->addF(total_force);
161  neighbour_particle->subF(total_force);
162  }
163  }
164  }
165  }
166  }
167 
168  // remove the periodic halo particles
169  particles.erase(particles.begin() + original_particle_length, particles.end());
171 }
std::vector< Cell * > & getNeighbourReferences()
Get the reference vector for the neighbouring cells.
Definition: Cell.cpp:9
std::vector< std::vector< Cell * > > iteration_orders
std::vector< Particle >::iterator end() override
The end iterator for the internal data structure.
Class to represent a particle.
Definition: Particle.h:26
bool isLocked() const
Gets whether the particle is locked in space.
Definition: Particle.h:202
void subF(const std::array< double, 3 > &force)
Subtracts a force from the particle.
Definition: Particle.h:140
void addF(const std::array< double, 3 > &force)
Adds a force to the particle.
Definition: Particle.h:129
static size_t applyBoundaryConditions(LinkedCellsContainer &container)
Applies the boundary conditions for the periodic boundary condition.
static void applyBoundaryConditions(LinkedCellsContainer &container)
Applies the boundary conditions for the reflective boundary condition.
auto L2NormSquared(const Container &c)
Definition: ArrayUtils.h:194

◆ applySimpleForces()

void LinkedCellsContainer::applySimpleForces ( const std::vector< std::shared_ptr< SimpleForceSource >> &  simple_force_sources)
overridevirtual

Applies the given simple force sources to the particles.

Parameters
simple_force_sourcesList of simple force sources to be applied

Applies the given simple force sources to the particles in the container.

Implements ParticleContainer.

Definition at line 108 of file LinkedCellsContainer.cpp.

108  {
109  for (Particle& p : particles) {
110  if (p.isLocked()) continue;
111  for (const auto& force_source : simple_force_sources) {
112  p.setF(p.getF() + force_source->calculateForce(p));
113  }
114  }
115 }

◆ applyTargettedForces()

void LinkedCellsContainer::applyTargettedForces ( const std::vector< std::shared_ptr< TargettedForceSource >> &  targetted_force_sources,
double  curr_simulation_time 
)
overridevirtual

Applies the given targetted force sources to the particles.

Parameters
targetted_force_sourcesList of targetted force sources to be applied
curr_simulation_timeCurrent simulation time

Applies the given targetted force sources to the particles in the container.

Implements ParticleContainer.

Definition at line 173 of file LinkedCellsContainer.cpp.

174  {
175  for (const auto& force_source : force_sources) {
176  force_source->applyForce(particles, curr_simulation_time);
177  }
178 }

◆ begin() [1/2]

std::vector< Particle >::const_iterator LinkedCellsContainer::begin ( ) const
overridevirtual

The begin const iterator for the internal data structure.

Returns
Const iterator to the first particle

Implements ParticleContainer.

Definition at line 199 of file LinkedCellsContainer.cpp.

199 { return particles.begin(); }

◆ begin() [2/2]

std::vector< Particle >::iterator LinkedCellsContainer::begin ( )
overridevirtual

The begin iterator for the internal data structure.

Returns
Iterator to the first particle

Implements ParticleContainer.

Definition at line 195 of file LinkedCellsContainer.cpp.

195 { return particles.begin(); }

◆ boundaryConditionToString()

std::string LinkedCellsContainer::boundaryConditionToString ( const BoundaryCondition bc)
static

Returns a string description of a boundary condition.

Definition at line 240 of file LinkedCellsContainer.cpp.

240  {
241  switch (bc) {
243  return "Outflow";
245  return "Reflective";
247  return "Periodic";
248  default:
249  return "Unknown";
250  }
251 };

◆ cellCoordToCellIndex()

int LinkedCellsContainer::cellCoordToCellIndex ( int  cx,
int  cy,
int  cz 
) const

Maps the cell coordinates to the corresponding index in the internal cell vector.

Parameters
cxx coordinate of the cell (valid range: -1 to domain_num_cells[0])
cyy coordinate of the cell (valid range: -1 to domain_num_cells[1])
czz coordinate of the cell (valid range: -1 to domain_num_cells[2])
Returns
index of the cell if it exists, -1 otherwise

Definition at line 217 of file LinkedCellsContainer.cpp.

217  {
218  if (cx < -1 || cx > domain_num_cells[0] || cy < -1 || cy > domain_num_cells[1] || cz < -1 || cz > domain_num_cells[2]) {
219  return -1;
220  }
221  return (cx + 1) * (domain_num_cells[1] + 2) * (domain_num_cells[2] + 2) + (cy + 1) * (domain_num_cells[2] + 2) + (cz + 1);
222 }

◆ deleteHaloParticles()

void LinkedCellsContainer::deleteHaloParticles ( )
private

Removes all particles in the halo cells from the particles vector. ATTENTION: A particle reference update must be triggered manually after this operation!

Definition at line 410 of file LinkedCellsContainer.cpp.

410  {
411  for (Cell* cell : halo_cell_references) {
412  for (Particle* p : cell->getParticleReferences()) {
413  particles.erase(particles.begin() + (p - &particles[0]));
414  }
415  }
416 }

◆ end() [1/2]

std::vector< Particle >::const_iterator LinkedCellsContainer::end ( ) const
overridevirtual

The end const iterator for the internal data structure.

Returns
Const end iterator for this container

Implements ParticleContainer.

Definition at line 201 of file LinkedCellsContainer.cpp.

201 { return particles.end(); }

◆ end() [2/2]

std::vector< Particle >::iterator LinkedCellsContainer::end ( )
overridevirtual

The end iterator for the internal data structure.

Returns
Iterator to the end of the container

Implements ParticleContainer.

Definition at line 197 of file LinkedCellsContainer.cpp.

197 { return particles.end(); }

◆ getBoundaryCells()

const std::vector< Cell * > & LinkedCellsContainer::getBoundaryCells ( ) const

Returns the pointers of the boundary cells.

Returns
Vector of pointers to the boundary cells

Returns the pointers of the boundary cells in a vector.

Definition at line 211 of file LinkedCellsContainer.cpp.

211 { return boundary_cell_references; }

◆ getCells()

const std::vector< Cell > & LinkedCellsContainer::getCells ( )

Returns the cells.

Returns
Cells

Returns the cells as a vector of Cell objects.

Definition at line 209 of file LinkedCellsContainer.cpp.

209 { return cells; }

◆ getCellSize()

const std::array< double, 3 > & LinkedCellsContainer::getCellSize ( ) const

Returns the cell size.

Returns
Cell size

Returns the cell size as a 3D array.

Definition at line 213 of file LinkedCellsContainer.cpp.

213 { return cell_size; }

◆ getCutoffRadius()

double LinkedCellsContainer::getCutoffRadius ( ) const

Returns the cutoff radius.

Returns
Cutoff radius

Returns the cutoff radius used for the force calculation.

Definition at line 207 of file LinkedCellsContainer.cpp.

207 { return cutoff_radius; }

◆ getDomainNumCells()

const std::array< int, 3 > & LinkedCellsContainer::getDomainNumCells ( ) const

Returns the number of cells in each dimension.

Returns
Number of cells in each dimension

Returns the number of cells in each dimension as a 3D array.

Definition at line 215 of file LinkedCellsContainer.cpp.

215 { return domain_num_cells; }

◆ getDomainSize()

const std::array< double, 3 > & LinkedCellsContainer::getDomainSize ( ) const

Returns the domain size.

Returns
Domain size

Returns the domain size as a 3D array.

Definition at line 205 of file LinkedCellsContainer.cpp.

205 { return domain_size; }

◆ getParticles()

const std::vector< Particle > & LinkedCellsContainer::getParticles ( ) const
overridevirtual

Returns a vector of all particles in the container.

Returns
Vector of all particles in the container

Implements ParticleContainer.

Definition at line 203 of file LinkedCellsContainer.cpp.

203 { return particles; }

◆ initCellNeighbourReferences()

void LinkedCellsContainer::initCellNeighbourReferences ( )
private

Sets the neighbour references for each cell in the cell vector.

Definition at line 316 of file LinkedCellsContainer.cpp.

316  {
317  // Loop through each cell according to their cell coordinates
318  for (int cx = -1; cx < domain_num_cells[0] + 1; ++cx) {
319  for (int cy = -1; cy < domain_num_cells[1] + 1; ++cy) {
320  for (int cz = -1; cz < domain_num_cells[2] + 1; ++cz) {
321  Cell& cell = cells.at(cellCoordToCellIndex(cx, cy, cz));
322 
323  // Loop through each of the current cells neighbour cells according to their cell coordinates
324  // except the current cell itself
325  for (int dx = -1; dx < 2; ++dx) {
326  for (int dy = -1; dy < 2; ++dy) {
327  for (int dz = -1; dz < 2; ++dz) {
328  if (dx == 0 && dy == 0 && dz == 0) continue;
329 
330  // Get the cell index of the current neighbour cell
331  int cell_index = cellCoordToCellIndex(cx + dx, cy + dy, cz + dz);
332 
333  // If the neighbour cell would be out of bounds, skip it
334  if (cell_index == -1) continue;
335 
336  // Add the neighbour to the current cells neighbour references
337  Cell& curr_neighbour = cells.at(cell_index);
338  cell.addNeighbourReference(&curr_neighbour);
339  }
340  }
341  }
342  }
343  }
344  }
345 }
void addNeighbourReference(Cell *c)
Adds a neighbour reference to the cell.
Definition: Cell.cpp:15
int cellCoordToCellIndex(int cx, int cy, int cz) const
Maps the cell coordinates to the corresponding index in the internal cell vector.

◆ initCells()

void LinkedCellsContainer::initCells ( )
private

Populates the cell vector and sets the cells types.

Definition at line 257 of file LinkedCellsContainer.cpp.

257  {
258  for (int cx = -1; cx < domain_num_cells[0] + 1; ++cx) {
259  for (int cy = -1; cy < domain_num_cells[1] + 1; ++cy) {
260  for (int cz = -1; cz < domain_num_cells[2] + 1; ++cz) {
261  if (cx < 0 || cx >= domain_num_cells[0] || cy < 0 || cy >= domain_num_cells[1] || cz < 0 || cz >= domain_num_cells[2]) {
262  cells.emplace_back(Cell::CellType::HALO);
263  halo_cell_references.push_back(&cells.back());
264 
265  if (cx == -1) {
266  left_halo_cell_references.push_back(&cells.back());
267  }
268  if (cx == domain_num_cells[0]) {
269  right_halo_cell_references.push_back(&cells.back());
270  }
271  if (cy == -1) {
272  bottom_halo_cell_references.push_back(&cells.back());
273  }
274  if (cy == domain_num_cells[1]) {
275  top_halo_cell_references.push_back(&cells.back());
276  }
277  if (cz == -1) {
278  back_halo_cell_references.push_back(&cells.back());
279  }
280  if (cz == domain_num_cells[2]) {
281  front_halo_cell_references.push_back(&cells.back());
282  }
283  } else if (cx == 0 || cx == domain_num_cells[0] - 1 || cy == 0 || cy == domain_num_cells[1] - 1 || cz == 0 ||
284  cz == domain_num_cells[2] - 1) {
285  cells.emplace_back(Cell::CellType::BOUNDARY);
286  boundary_cell_references.push_back(&cells.back());
287  domain_cell_references.push_back(&cells.back());
288 
289  if (cx == 0) {
290  left_boundary_cell_references.push_back(&cells.back());
291  }
292  if (cx == domain_num_cells[0] - 1) {
293  right_boundary_cell_references.push_back(&cells.back());
294  }
295  if (cy == 0) {
296  bottom_boundary_cell_references.push_back(&cells.back());
297  }
298  if (cy == domain_num_cells[1] - 1) {
299  top_boundary_cell_references.push_back(&cells.back());
300  }
301  if (cz == 0) {
302  back_boundary_cell_references.push_back(&cells.back());
303  }
304  if (cz == domain_num_cells[2] - 1) {
305  front_boundary_cell_references.push_back(&cells.back());
306  }
307  } else {
308  cells.emplace_back(Cell::CellType::INNER);
309  domain_cell_references.push_back(&cells.back());
310  }
311  }
312  }
313  }
314 }
std::vector< Cell * > right_halo_cell_references
References to the halo cells on the right (x = domain_num_cells[0])
std::vector< Cell * > bottom_boundary_cell_references
References to the boundary cells on the bottom (y = 0)
std::vector< Cell * > top_boundary_cell_references
References to the boundary cells on the top (y = domain_num_cells[1]-1)
std::vector< Cell * > front_halo_cell_references
References to the halo cells on the front (z = domain_num_cells[2])
std::vector< Cell * > bottom_halo_cell_references
References to the halo cells on the bottom (y = -1)
std::vector< Cell * > left_halo_cell_references
References to the halo cells on the left (x = -1)
std::vector< Cell * > right_boundary_cell_references
References to the boundary cells on the right (x = domain_num_cells[0]-1)
std::vector< Cell * > top_halo_cell_references
References to the halo cells on the top (y = domain_num_cells[1])
std::vector< Cell * > left_boundary_cell_references
References to the boundary cells on the left (x = 0)
std::vector< Cell * > back_halo_cell_references
References to the halo cells on the back (z = -1)
std::vector< Cell * > back_boundary_cell_references
References to the boundary cells on the back (z = 0)
std::vector< Cell * > front_boundary_cell_references
References to the boundary cells on the front (z = domain_num_cells[2]-1)

◆ initIterationOrders()

void LinkedCellsContainer::initIterationOrders ( )
private

Sets the iteration orders for parallel execution.

Definition at line 347 of file LinkedCellsContainer.cpp.

347  {
348 #if PARALLEL_V2
349  Logger::logger->warn("Creating iteration orders for parallel v2");
350  std::vector<Cell*> iteration_order;
351 
352  for (size_t num = 0; num < cells.size(); ++num) {
353  iteration_order.push_back(&cells[num]);
354  }
355 
356  // shuffle the iteration order
357  std::random_device rd;
358  std::mt19937 g(rd());
359  std::shuffle(iteration_order.begin(), iteration_order.end(), g);
360 
361  iteration_orders.push_back(iteration_order);
362 #else
363  Logger::logger->warn("Creating iteration orders for parallel v1 (c_18 traversal)");
364 
365  std::vector<std::array<int, 3>> start_offsets;
366 
367  int d_x = 2;
368  int d_y = 3;
369  int d_z = 3;
370 
371  for (int x = 0; x < d_x; ++x) {
372  for (int y = 0; y < d_y; ++y) {
373  for (int z = 0; z < d_z; ++z) {
374  start_offsets.push_back({x - 1, y - 1, z - 1});
375  }
376  }
377  }
378 
379  for (size_t i = 0; i < start_offsets.size(); ++i) {
380  std::vector<Cell*> iteration_order;
381  const std::array<int, 3>& start_offset = start_offsets[i];
382 
383  for (int cx = start_offset[0]; cx <= domain_num_cells[0]; cx += d_x) {
384  for (int cy = start_offset[1]; cy <= domain_num_cells[1]; cy += d_y) {
385  for (int cz = start_offset[2]; cz <= domain_num_cells[2]; cz += d_z) {
386  iteration_order.push_back(&cells.at(cellCoordToCellIndex(cx, cy, cz)));
387  }
388  }
389  }
390 
391  iteration_orders.push_back(iteration_order);
392  }
393 #endif
394 }

◆ operator[]()

Particle & LinkedCellsContainer::operator[] ( int  i)
overridevirtual

Overload of the [] operator to access the particles in the container.

Parameters
iIndex of the particle
Returns
Particle

Implements ParticleContainer.

Definition at line 193 of file LinkedCellsContainer.cpp.

193 { return particles[i]; }

◆ particlePosToCell() [1/2]

Cell * LinkedCellsContainer::particlePosToCell ( const std::array< double, 3 > &  pos)

Maps the particle position to the corresponding cell index in the internal cell vector.

Parameters
posPosition of the particle
Returns
index of the cell if it exists, -1 otherwise

Definition at line 224 of file LinkedCellsContainer.cpp.

224 { return particlePosToCell(pos[0], pos[1], pos[2]); }

◆ particlePosToCell() [2/2]

Cell * LinkedCellsContainer::particlePosToCell ( double  x,
double  y,
double  z 
)

Maps the particle position to the corresponding cell index in the internal cell vector.

Parameters
xx coordinate of the particle
yy coordinate of the particle
zz coordinate of the particle
Returns
index of the cell if it exists, -1 otherwise

Definition at line 226 of file LinkedCellsContainer.cpp.

226  {
227  int cx = static_cast<int>(std::floor(x / cell_size[0]));
228  int cy = static_cast<int>(std::floor(y / cell_size[1]));
229  int cz = static_cast<int>(std::floor(z / cell_size[2]));
230 
231  int cell_index = cellCoordToCellIndex(cx, cy, cz);
232  if (cell_index == -1) {
233  Logger::logger->error("Particle is outside of cells. Position: [{}, {}, {}]", x, y, z);
234  throw std::runtime_error("A particle is outside of the cells");
235  }
236 
237  return &cells[cell_index];
238 }

◆ prepareForceCalculation()

void LinkedCellsContainer::prepareForceCalculation ( )
overridevirtual

Prepares everything for the force calculations (must be called before applySimpleForces and applyPairwiseForces)

Implements ParticleContainer.

Definition at line 94 of file LinkedCellsContainer.cpp.

94  {
95  // update the particle references in the cells in case the particles have moved
97 
98  // delete particles from halo cells
100 
101  // teleport particles to other side of container
103 
104  // update the particle references in the cells in case the particles have moved
106 }
static void pre(LinkedCellsContainer &container)
Applies the preconditioning step for the outflow boundary condition.
static void pre(LinkedCellsContainer &container)
Applies the preconditioning step for the periodic boundary condition.

◆ reserve()

void LinkedCellsContainer::reserve ( size_t  n)
overridevirtual

Reserves space for n particles. This is useful if the number of particles is known in advance and prevents reallocation of memory for the internal dynamic array of particles, when inserting new particles.

Parameters
nTotal amount of particles expected to store in the container

Implements ParticleContainer.

Definition at line 180 of file LinkedCellsContainer.cpp.

180  {
181  Logger::logger->debug("Reserving space for {} particles", n);
182 
183  size_t old_capacity = particles.capacity();
184  particles.reserve(n);
185 
186  if (old_capacity != particles.capacity()) {
188  }
189 }

◆ size()

size_t LinkedCellsContainer::size ( ) const
overridevirtual

Returns the number of particles in the container.

Returns
Number of particles in the container

Implements ParticleContainer.

Definition at line 191 of file LinkedCellsContainer.cpp.

191 { return particles.size(); }

◆ updateCellsParticleReferences()

void LinkedCellsContainer::updateCellsParticleReferences ( )
private

Updates the particle references in the cells. This is necessary after a reallocation of the internal particle vector.

Definition at line 396 of file LinkedCellsContainer.cpp.

396  {
397  // clear the particle references in the cells
398  for (Cell& cell : cells) {
399  cell.clearParticleReferences();
400  }
401 
402  // add the particle references to the cells
403  for (Particle& p : particles) {
404  Cell* cell = particlePosToCell(p.getX());
405 
406  cell->addParticleReference(&p);
407  }
408 }

Friends And Related Function Documentation

◆ OutflowBoundaryType

friend class OutflowBoundaryType
friend

Friend class to allow access to the internal data structures.

Definition at line 287 of file LinkedCellsContainer.h.

◆ PeriodicBoundaryType

friend class PeriodicBoundaryType
friend

Friend class to allow access to the internal data structures.

Definition at line 292 of file LinkedCellsContainer.h.

◆ ReflectiveBoundaryType

friend class ReflectiveBoundaryType
friend

Friend class to allow access to the internal data structures.

Definition at line 282 of file LinkedCellsContainer.h.

Member Data Documentation

◆ back_boundary_cell_references

std::vector<Cell*> LinkedCellsContainer::back_boundary_cell_references
private

References to the boundary cells on the back (z = 0)

Definition at line 370 of file LinkedCellsContainer.h.

◆ back_halo_cell_references

std::vector<Cell*> LinkedCellsContainer::back_halo_cell_references
private

References to the halo cells on the back (z = -1)

Definition at line 402 of file LinkedCellsContainer.h.

◆ bottom_boundary_cell_references

std::vector<Cell*> LinkedCellsContainer::bottom_boundary_cell_references
private

References to the boundary cells on the bottom (y = 0)

Definition at line 360 of file LinkedCellsContainer.h.

◆ bottom_halo_cell_references

std::vector<Cell*> LinkedCellsContainer::bottom_halo_cell_references
private

References to the halo cells on the bottom (y = -1)

Definition at line 392 of file LinkedCellsContainer.h.

◆ boundary_cell_references

std::vector<Cell*> LinkedCellsContainer::boundary_cell_references
private

References to the boundary cells.

Definition at line 338 of file LinkedCellsContainer.h.

◆ boundary_types

std::array<BoundaryCondition, 6> LinkedCellsContainer::boundary_types
private

The boundary types for each side of the domain (order in array: left, right, bottom, top, back, front)

Definition at line 313 of file LinkedCellsContainer.h.

◆ cell_size

std::array<double, 3> LinkedCellsContainer::cell_size
private

Cell size in each dimension.

Definition at line 318 of file LinkedCellsContainer.h.

◆ cells

std::vector<Cell> LinkedCellsContainer::cells
private

Internal data structure for the cells.

Definition at line 328 of file LinkedCellsContainer.h.

◆ cutoff_radius

double LinkedCellsContainer::cutoff_radius
private

Cutoff radius for the force calculation.

Definition at line 308 of file LinkedCellsContainer.h.

◆ domain_cell_references

std::vector<Cell*> LinkedCellsContainer::domain_cell_references
private

References to the domain cells.

Definition at line 333 of file LinkedCellsContainer.h.

◆ domain_num_cells

std::array<int, 3> LinkedCellsContainer::domain_num_cells
private

Number of cells in each dimension.

Definition at line 323 of file LinkedCellsContainer.h.

◆ domain_size

std::array<double, 3> LinkedCellsContainer::domain_size
private

Domain size in each dimension.

Definition at line 303 of file LinkedCellsContainer.h.

◆ front_boundary_cell_references

std::vector<Cell*> LinkedCellsContainer::front_boundary_cell_references
private

References to the boundary cells on the front (z = domain_num_cells[2]-1)

Definition at line 375 of file LinkedCellsContainer.h.

◆ front_halo_cell_references

std::vector<Cell*> LinkedCellsContainer::front_halo_cell_references
private

References to the halo cells on the front (z = domain_num_cells[2])

Definition at line 407 of file LinkedCellsContainer.h.

◆ halo_cell_references

std::vector<Cell*> LinkedCellsContainer::halo_cell_references
private

References to the halo cells.

Definition at line 343 of file LinkedCellsContainer.h.

◆ iteration_orders

std::vector<std::vector<Cell*> > LinkedCellsContainer::iteration_orders
private

Definition at line 15 of file LinkedCellsContainer.h.

◆ left_boundary_cell_references

std::vector<Cell*> LinkedCellsContainer::left_boundary_cell_references
private

References to the boundary cells on the left (x = 0)

Definition at line 350 of file LinkedCellsContainer.h.

◆ left_halo_cell_references

std::vector<Cell*> LinkedCellsContainer::left_halo_cell_references
private

References to the halo cells on the left (x = -1)

Definition at line 382 of file LinkedCellsContainer.h.

◆ particles

std::vector<Particle> LinkedCellsContainer::particles
private

Internal data structure for the particles.

Definition at line 298 of file LinkedCellsContainer.h.

◆ right_boundary_cell_references

std::vector<Cell*> LinkedCellsContainer::right_boundary_cell_references
private

References to the boundary cells on the right (x = domain_num_cells[0]-1)

Definition at line 355 of file LinkedCellsContainer.h.

◆ right_halo_cell_references

std::vector<Cell*> LinkedCellsContainer::right_halo_cell_references
private

References to the halo cells on the right (x = domain_num_cells[0])

Definition at line 387 of file LinkedCellsContainer.h.

◆ top_boundary_cell_references

std::vector<Cell*> LinkedCellsContainer::top_boundary_cell_references
private

References to the boundary cells on the top (y = domain_num_cells[1]-1)

Definition at line 365 of file LinkedCellsContainer.h.

◆ top_halo_cell_references

std::vector<Cell*> LinkedCellsContainer::top_halo_cell_references
private

References to the halo cells on the top (y = domain_num_cells[1])

Definition at line 397 of file LinkedCellsContainer.h.


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