ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
Loading...
Searching...
No Matches
CellStructure Class Reference

Describes a cell structure / cell system. More...

#include <CellStructure.hpp>

+ Inheritance diagram for CellStructure:
+ Collaboration diagram for CellStructure:

Classes

struct  AoSoA_pack
 

Public Types

using ForceType = Kokkos::View< double **[3], Kokkos::LayoutRight >
 
using VirialType = Kokkos::View< double *[3], Kokkos::LayoutRight >
 
using memory_space = Kokkos::HostSpace
 
using ListAlgorithm = Cabana::HalfNeighborTag
 
using ListType = CustomVerletList< Kokkos::HostSpace, ListAlgorithm, Cabana::VerletLayout2D, Cabana::TeamVectorOpTag >
 

Public Member Functions

 CellStructure (BoxGeometry const &box)
 
virtual ~CellStructure ()
 
void update_particle_index (int id, Particle *p)
 Update local particle index.
 
void update_particle_index (Particle &p)
 Update local particle index.
 
void update_particle_index (ParticleList &pl)
 Update local particle index.
 
void clear_particle_index ()
 Clear the particles index.
 
Particleget_local_particle (int id)
 Get a local particle by id.
 
const Particleget_local_particle (int id) const
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
template<class InputRange , class OutputIterator >
void get_local_particles (InputRange ids, OutputIterator out)
 
CellStructureType decomposition_type () const
 
Utils::Vector3d max_cutoff () const
 Maximal cutoff supported by current cell system.
 
Utils::Vector3d max_range () const
 Maximal pair range supported by current cell system.
 
ParticleRange local_particles () const
 
ParticleRange ghost_particles () const
 
std::size_t count_local_particles () const
 
bool use_parallel_for_each_local_particle () const
 whether to use parallel version of for_each_local_particle
 
void for_each_local_particle (ParticleUnaryOp &&f, bool parallel=true) const
 Run a kernel on all local particles.
 
void for_each_ghost_particle (ParticleUnaryOp &&f) const
 Run a kernel on all ghost particles.
 
Particleadd_particle (Particle &&p)
 Add a particle.
 
Particleadd_local_particle (Particle &&p)
 Add a particle.
 
void remove_particle (int id)
 Remove a particle.
 
int get_max_local_particle_id () const
 Get the maximal particle id on this node.
 
int get_cached_max_local_particle_id () const
 
void remove_all_particles ()
 Remove all particles from the cell system.
 
ParticleDecomposition const & decomposition () const
 Get the underlying particle decomposition.
 
void set_resort_particles (Cells::Resort level)
 Increase the local resort level at least to level.
 
unsigned get_resort_particles () const
 Get the currently scheduled resort level.
 
void clear_resort_particles ()
 Set the resort level to sorted.
 
bool check_resort_required (Utils::Vector3d const &additional_offset={}) const
 Check whether a particle has moved further than half the skin since the last Verlet list update, thus requiring a resort.
 
auto get_le_pos_offset_at_last_resort () const
 
void ghosts_count ()
 Synchronize number of ghosts.
 
void ghosts_update (unsigned data_parts)
 Update ghost particles.
 
void update_ghosts_and_resort_particle (unsigned data_parts)
 Update ghost particles, with particle resort if needed.
 
void ghosts_reduce_forces ()
 Add forces from ghost particles to real particles.
 
void ghosts_reduce_rattle_correction ()
 Add rattle corrections from ghost particles to real particles.
 
void resort_particles (bool global_flag)
 Resort particles.
 
auto is_verlet_skin_set () const
 Whether the Verlet skin is set.
 
auto get_verlet_skin () const
 Get the Verlet skin.
 
void set_verlet_skin (double value)
 Set the Verlet skin.
 
void set_verlet_skin_heuristic ()
 Set the Verlet skin using a heuristic.
 
void update_verlet_stats (int n_steps, int n_verlet_updates)
 
auto get_verlet_reuse () const
 Average number of integration steps the Verlet list was re-used.
 
auto resolve_bond_partners (std::span< const int > partner_ids)
 Resolve ids to particles.
 
void set_atom_decomposition ()
 Set the particle decomposition to AtomDecomposition.
 
void set_regular_decomposition (double range, std::optional< std::pair< int, int > > fully_connected_boundary)
 Set the particle decomposition to RegularDecomposition.
 
void set_hybrid_decomposition (double cutoff_regular, std::set< int > n_square_types)
 Set the particle decomposition to HybridDecomposition.
 
auto get_max_id () const
 
void set_kokkos_handle (std::shared_ptr< KokkosHandle > handle)
 
void rebuild_local_properties (double pair_cutoff)
 
void reset_local_properties ()
 
void reset_local_force ()
 
auto & get_id_to_index ()
 
auto & get_local_force ()
 
auto & get_local_torque ()
 
auto & get_local_virial ()
 
auto & get_aosoa ()
 
auto const & get_unique_particles () const
 
auto const & get_verlet_list_cabana () const
 
void clear_local_properties ()
 
auto is_verlet_list_cabana_rebuild_needed () const
 
auto prepare_verlet_list_cabana (double cutoff)
 Reset local properties of the Verlet list.
 
void rebuild_verlet_list_cabana (auto &&kernel, bool rebuild_verlet_list)
 
void set_index_map ()
 
void cell_list_loop (auto &&kernel)
 
template<class BondKernel >
void bond_loop (BondKernel const &bond_kernel)
 Bonded pair loop.
 
template<class PairKernel >
void non_bonded_loop (PairKernel pair_kernel)
 Non-bonded pair loop.
 
template<class PairKernel , class VerletCriterion >
void non_bonded_loop (PairKernel pair_kernel, const VerletCriterion &verlet_criterion)
 Non-bonded pair loop with potential use of verlet lists.
 
void check_particle_index () const
 Check that particle index is commensurate with particles.
 
void check_particle_sorting () const
 Check that particles are in the correct cell.
 
Cellfind_current_cell (const Particle &p)
 Find cell a particle is stored in.
 
template<class Kernel >
bool run_on_particle_short_range_neighbors (Particle const &p, Kernel &kernel)
 Run kernel on all particles inside local cell and its neighbors.
 
- Public Member Functions inherited from System::Leaf< CellStructure >
void bind_system (std::shared_ptr< System > const &system)
 
void detach_system (std::shared_ptr< System > const &system)
 

Public Attributes

bool use_verlet_list = true
 

Static Public Attributes

static constexpr auto vector_length = 1
 

Additional Inherited Members

- Protected Member Functions inherited from System::Leaf< CellStructure >
auto & get_system ()
 
auto & get_system () const
 
- Protected Attributes inherited from System::Leaf< CellStructure >
std::weak_ptr< Systemm_system
 

Detailed Description

Describes a cell structure / cell system.

Contains information about the communication of cell contents (particles, ghosts, ...) between different nodes and the relation between particle positions and the cell system. All other properties of the cell system which are not common between different cell systems have to be stored in separate structures.

Definition at line 169 of file CellStructure.hpp.

Member Typedef Documentation

◆ ForceType

using CellStructure::ForceType = Kokkos::View<double **[3], Kokkos::LayoutRight>

Definition at line 174 of file CellStructure.hpp.

◆ ListAlgorithm

using CellStructure::ListAlgorithm = Cabana::HalfNeighborTag

Definition at line 177 of file CellStructure.hpp.

◆ ListType

using CellStructure::ListType = CustomVerletList<Kokkos::HostSpace, ListAlgorithm, Cabana::VerletLayout2D, Cabana::TeamVectorOpTag>

Definition at line 178 of file CellStructure.hpp.

◆ memory_space

using CellStructure::memory_space = Kokkos::HostSpace

Definition at line 176 of file CellStructure.hpp.

◆ VirialType

using CellStructure::VirialType = Kokkos::View<double *[3], Kokkos::LayoutRight>

Definition at line 175 of file CellStructure.hpp.

Constructor & Destructor Documentation

◆ CellStructure()

CellStructure::CellStructure ( BoxGeometry const &  box)

Definition at line 227 of file CellStructure.cpp.

◆ ~CellStructure()

CellStructure::~CellStructure ( )
virtual

Definition at line 75 of file CellStructure.cpp.

References clear_local_properties().

Member Function Documentation

◆ add_local_particle()

Particle * CellStructure::add_local_particle ( Particle &&  p)

Add a particle.

Moves a particle into the cell system, if it belongs to this node. Otherwise this does not have an effect and the particle is discarded. This can be used to add a particle without knowledge where it should be placed by calling the function on all nodes, it will then add the particle in exactly one place.

Parameters
pParticle to add.
Returns
Pointer to particle if it is local, null otherwise.

Definition at line 300 of file CellStructure.cpp.

◆ add_particle()

Particle * CellStructure::add_particle ( Particle &&  p)

Add a particle.

Moves a particle into the cell system. This adds a particle to the local node, irrespective of where it belongs.

Parameters
pParticle to add.
Returns
Pointer to the particle in the cell system.

Definition at line 310 of file CellStructure.cpp.

References decomposition(), ParticleDecomposition::local_cells(), Cells::RESORT_GLOBAL, Cells::RESORT_LOCAL, and set_resort_particles().

Referenced by Mpiio::mpi_mpiio_common_read(), and CollisionDetection::place_vs_and_relate_to_particle().

◆ bond_loop()

template<class BondKernel >
void CellStructure::bond_loop ( BondKernel const &  bond_kernel)
inline

Bonded pair loop.

Parameters
bond_kernelKernel to apply

Definition at line 828 of file CellStructure.hpp.

References local_particles().

Referenced by add_oif_global_forces(), cabana_short_range(), calc_oif_mesh(), and compute_correction_vector().

◆ cell_list_loop()

void CellStructure::cell_list_loop ( auto &&  kernel)
inline

Definition at line 777 of file CellStructure.hpp.

Referenced by cabana_short_range().

◆ check_particle_index()

void CellStructure::check_particle_index ( ) const

Check that particle index is commensurate with particles.

For each local particles is checked that has a correct entry in the particles index, and that there are no excess (non-existing) particles in the index.

Definition at line 230 of file CellStructure.cpp.

References get_local_particle(), get_max_local_particle_id(), Particle::id(), and local_particles().

Referenced by resort_particles().

◆ check_particle_sorting()

void CellStructure::check_particle_sorting ( ) const

Check that particles are in the correct cell.

This checks for all local particles that the result of particles_to_cell is the cell the particles is actually in, e.g. that the particles are sorted according to particles_to_cell.

Definition at line 263 of file CellStructure.cpp.

References decomposition(), and ParticleDecomposition::local_cells().

Referenced by resort_particles().

◆ check_resort_required()

bool CellStructure::check_resort_required ( Utils::Vector3d const &  additional_offset = {}) const

Check whether a particle has moved further than half the skin since the last Verlet list update, thus requiring a resort.

Parameters
additional_offsetOffset which is added to the distance the particle has travelled when comparing to half the Verlet skin (e.g., for Lees-Edwards BC).
Returns
Whether a resort is needed.

Definition at line 522 of file CellStructure.cpp.

References Utils::Vector< T, N >::norm2(), reduce_over_local_particles(), and Utils::sqr().

Referenced by correct_position_shake(), and vs_relative_update_particles().

◆ clear_local_properties()

void CellStructure::clear_local_properties ( )

Definition at line 84 of file CellStructure.cpp.

Referenced by ~CellStructure().

◆ clear_particle_index()

void CellStructure::clear_particle_index ( )
inline

Clear the particles index.

Definition at line 272 of file CellStructure.hpp.

◆ clear_resort_particles()

void CellStructure::clear_resort_particles ( )
inline

Set the resort level to sorted.

Definition at line 505 of file CellStructure.hpp.

References Cells::RESORT_NONE.

Referenced by update_ghosts_and_resort_particle().

◆ count_local_particles()

std::size_t CellStructure::count_local_particles ( ) const
inline

Definition at line 346 of file CellStructure.hpp.

Referenced by set_index_map().

◆ decomposition()

ParticleDecomposition const & CellStructure::decomposition ( ) const
inline

◆ decomposition_type()

CellStructureType CellStructure::decomposition_type ( ) const
inline

Definition at line 330 of file CellStructure.hpp.

◆ find_current_cell()

Cell * CellStructure::find_current_cell ( const Particle p)
inline

Find cell a particle is stored in.

For local particles, this returns the cell they are stored in, otherwise nullptr is returned.

Parameters
pParticle to find cell for
Returns
Cell for particle or nullptr.

Definition at line 886 of file CellStructure.hpp.

References get_resort_particles(), and Particle::is_ghost().

Referenced by run_on_particle_short_range_neighbors().

◆ for_each_ghost_particle()

void CellStructure::for_each_ghost_particle ( ParticleUnaryOp &&  f) const
inline

Run a kernel on all ghost particles.

The kernel is assumed to be thread-safe.

Definition at line 384 of file CellStructure.hpp.

References ghost_particles().

Referenced by force_calc_icc(), and init_forces_ghosts().

◆ for_each_local_particle()

void CellStructure::for_each_local_particle ( ParticleUnaryOp &&  f,
bool  parallel = true 
) const
inline

◆ get_aosoa()

auto & CellStructure::get_aosoa ( )
inline

◆ get_cached_max_local_particle_id()

int CellStructure::get_cached_max_local_particle_id ( ) const
inline

◆ get_id_to_index()

auto & CellStructure::get_id_to_index ( )
inline

◆ get_le_pos_offset_at_last_resort()

auto CellStructure::get_le_pos_offset_at_last_resort ( ) const
inline

Definition at line 518 of file CellStructure.hpp.

◆ get_local_force()

auto & CellStructure::get_local_force ( )
inline

◆ get_local_particle() [1/2]

◆ get_local_particle() [2/2]

const Particle * CellStructure::get_local_particle ( int  id) const
inline

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

Definition at line 315 of file CellStructure.hpp.

◆ get_local_particles()

template<class InputRange , class OutputIterator >
void CellStructure::get_local_particles ( InputRange  ids,
OutputIterator  out 
)
inline

Definition at line 325 of file CellStructure.hpp.

References get_local_particle().

Referenced by resolve_bond_partners().

◆ get_local_torque()

auto & CellStructure::get_local_torque ( )
inline

Definition at line 733 of file CellStructure.hpp.

Referenced by rebuild_local_properties(), and reset_local_properties().

◆ get_local_virial()

auto & CellStructure::get_local_virial ( )
inline

Definition at line 736 of file CellStructure.hpp.

Referenced by reset_local_properties().

◆ get_max_id()

auto CellStructure::get_max_id ( ) const
inline

Definition at line 723 of file CellStructure.hpp.

◆ get_max_local_particle_id()

int CellStructure::get_max_local_particle_id ( ) const

Get the maximal particle id on this node.

This returns the highest particle id on this node, or -1 if there are no particles on this node.

Definition at line 324 of file CellStructure.cpp.

Referenced by check_particle_index().

◆ get_resort_particles()

unsigned CellStructure::get_resort_particles ( ) const
inline

Get the currently scheduled resort level.

Definition at line 500 of file CellStructure.hpp.

Referenced by cabana_short_range(), and find_current_cell().

◆ get_unique_particles()

auto const & CellStructure::get_unique_particles ( ) const
inline

◆ get_verlet_list_cabana()

auto const & CellStructure::get_verlet_list_cabana ( ) const
inline

Definition at line 740 of file CellStructure.hpp.

Referenced by cabana_short_range().

◆ get_verlet_reuse()

auto CellStructure::get_verlet_reuse ( ) const
inline

Average number of integration steps the Verlet list was re-used.

Definition at line 586 of file CellStructure.hpp.

◆ get_verlet_skin()

auto CellStructure::get_verlet_skin ( ) const
inline

Get the Verlet skin.

Definition at line 569 of file CellStructure.hpp.

Referenced by lb_tracers_propagate().

◆ ghost_particles()

◆ ghosts_count()

void CellStructure::ghosts_count ( )

◆ ghosts_reduce_forces()

void CellStructure::ghosts_reduce_forces ( )

◆ ghosts_reduce_rattle_correction()

void CellStructure::ghosts_reduce_rattle_correction ( )

Add rattle corrections from ghost particles to real particles.

Definition at line 370 of file CellStructure.cpp.

References decomposition(), System::Leaf< CellStructure >::get_system(), ghost_communicator(), and GHOSTTRANS_RATTLE.

Referenced by correct_position_shake(), and correct_velocity_shake().

◆ ghosts_update()

void CellStructure::ghosts_update ( unsigned  data_parts)

Update ghost particles.

Update ghost particles with data from the real particles.

Parameters
data_partsParticle parts to update, combination of Cells::DataPart

Definition at line 361 of file CellStructure.cpp.

References decomposition(), System::Leaf< CellStructure >::get_system(), ghost_communicator(), and map_data_parts().

Referenced by correct_position_shake(), correct_velocity_shake(), update_ghosts_and_resort_particle(), vs_com_update_particles(), and vs_relative_update_particles().

◆ is_verlet_list_cabana_rebuild_needed()

auto CellStructure::is_verlet_list_cabana_rebuild_needed ( ) const
inline

Definition at line 743 of file CellStructure.hpp.

Referenced by prepare_verlet_list_cabana(), and rebuild_verlet_list_cabana().

◆ is_verlet_skin_set()

auto CellStructure::is_verlet_skin_set ( ) const
inline

Whether the Verlet skin is set.

Definition at line 566 of file CellStructure.hpp.

Referenced by set_verlet_skin_heuristic().

◆ local_particles()

◆ max_cutoff()

Utils::Vector3d CellStructure::max_cutoff ( ) const
inline

Maximal cutoff supported by current cell system.

Definition at line 333 of file CellStructure.hpp.

References decomposition(), and ParticleDecomposition::max_cutoff().

Referenced by set_verlet_skin_heuristic().

◆ max_range()

Utils::Vector3d CellStructure::max_range ( ) const
inline

Maximal pair range supported by current cell system.

Definition at line 336 of file CellStructure.hpp.

References decomposition(), and ParticleDecomposition::max_range().

Referenced by set_verlet_skin_heuristic().

◆ non_bonded_loop() [1/2]

template<class PairKernel >
void CellStructure::non_bonded_loop ( PairKernel  pair_kernel)
inline

Non-bonded pair loop.

Parameters
pair_kernelKernel to apply

Definition at line 837 of file CellStructure.hpp.

Referenced by force_calc_icc().

◆ non_bonded_loop() [2/2]

template<class PairKernel , class VerletCriterion >
void CellStructure::non_bonded_loop ( PairKernel  pair_kernel,
const VerletCriterion verlet_criterion 
)
inline

Non-bonded pair loop with potential use of verlet lists.

Parameters
pair_kernelKernel to apply
verlet_criterionFilter for verlet lists.

Definition at line 847 of file CellStructure.hpp.

References use_verlet_list.

◆ prepare_verlet_list_cabana()

auto CellStructure::prepare_verlet_list_cabana ( double  cutoff)
inline

Reset local properties of the Verlet list.

Parameters
cutoffPair interaction cutoff.
Returns
True if a rebuild is needed.

Definition at line 752 of file CellStructure.hpp.

References is_verlet_list_cabana_rebuild_needed(), rebuild_local_properties(), reset_local_properties(), and set_index_map().

Referenced by update_cabana_state().

◆ rebuild_local_properties()

◆ rebuild_verlet_list_cabana()

void CellStructure::rebuild_verlet_list_cabana ( auto &&  kernel,
bool  rebuild_verlet_list 
)
inline

Definition at line 766 of file CellStructure.hpp.

References is_verlet_list_cabana_rebuild_needed().

Referenced by update_cabana_state().

◆ remove_all_particles()

void CellStructure::remove_all_particles ( )

Remove all particles from the cell system.

This allows linear time removal of all particles from the system, removing each particle individually would be quadratic.

Definition at line 331 of file CellStructure.cpp.

References decomposition(), and ParticleDecomposition::local_cells().

Referenced by Mpiio::mpi_mpiio_common_read().

◆ remove_particle()

void CellStructure::remove_particle ( int  id)

Remove a particle.

Removes a particle and all bonds pointing to it. This is a collective call.

Parameters
idId of particle to remove.

Definition at line 274 of file CellStructure.cpp.

References Utils::contains(), decomposition(), ParticleDecomposition::local_cells(), and update_particle_index().

◆ reset_local_force()

void CellStructure::reset_local_force ( )

Definition at line 170 of file CellStructure.cpp.

References get_local_force().

Referenced by force_calc_icc().

◆ reset_local_properties()

void CellStructure::reset_local_properties ( )

◆ resolve_bond_partners()

auto CellStructure::resolve_bond_partners ( std::span< const int >  partner_ids)
inline

Resolve ids to particles.

Exceptions
BondResolutionErrorif one of the ids was not found.
Parameters
partner_idsIds to resolve.
Returns
Vector of Particle pointers.

Definition at line 597 of file CellStructure.hpp.

References get_local_particles().

◆ resort_particles()

void CellStructure::resort_particles ( bool  global_flag)

◆ run_on_particle_short_range_neighbors()

template<class Kernel >
bool CellStructure::run_on_particle_short_range_neighbors ( Particle const &  p,
Kernel &  kernel 
)
inline

Run kernel on all particles inside local cell and its neighbors.

Parameters
pParticle to find cell for
kernelFunction with signature double(Particle const&, Particle const&, Utils::Vector3d const&)
Returns
false if cell is not found, otherwise true

Definition at line 905 of file CellStructure.hpp.

References ParticleDecomposition::box(), decomposition(), find_current_cell(), and ParticleDecomposition::minimum_image_distance().

◆ set_atom_decomposition()

void CellStructure::set_atom_decomposition ( )

Set the particle decomposition to AtomDecomposition.

Definition at line 413 of file CellStructure.cpp.

References comm_cart, System::Leaf< CellStructure >::get_system(), and NSQUARE.

◆ set_hybrid_decomposition()

void CellStructure::set_hybrid_decomposition ( double  cutoff_regular,
std::set< int >  n_square_types 
)

Set the particle decomposition to HybridDecomposition.

Parameters
cutoff_regularInteraction cutoff_regular.
n_square_typesParticle types to put into n_square decomposition.

Definition at line 436 of file CellStructure.cpp.

References comm_cart, System::Leaf< CellStructure >::get_system(), and HYBRID.

◆ set_index_map()

void CellStructure::set_index_map ( )

◆ set_kokkos_handle()

void CellStructure::set_kokkos_handle ( std::shared_ptr< KokkosHandle handle)

Definition at line 98 of file CellStructure.cpp.

◆ set_regular_decomposition()

void CellStructure::set_regular_decomposition ( double  range,
std::optional< std::pair< int, int > >  fully_connected_boundary 
)

Set the particle decomposition to RegularDecomposition.

Parameters
rangeInteraction range.
fully_connected_boundaryneighbor cell directions for Lees-Edwards.

Definition at line 424 of file CellStructure.cpp.

References comm_cart, System::Leaf< CellStructure >::get_system(), and REGULAR.

◆ set_resort_particles()

void CellStructure::set_resort_particles ( Cells::Resort  level)
inline

Increase the local resort level at least to level.

Definition at line 492 of file CellStructure.hpp.

Referenced by add_particle(), correct_position_shake(), lb_tracers_propagate(), ScriptInterface::System::rotate_system(), and vs_relative_update_particles().

◆ set_verlet_skin()

void CellStructure::set_verlet_skin ( double  value)

Set the Verlet skin.

Definition at line 450 of file CellStructure.cpp.

References System::Leaf< CellStructure >::get_system().

Referenced by set_verlet_skin_heuristic().

◆ set_verlet_skin_heuristic()

void CellStructure::set_verlet_skin_heuristic ( )

Set the Verlet skin using a heuristic.

Definition at line 458 of file CellStructure.cpp.

References System::Leaf< CellStructure >::get_system(), is_verlet_skin_set(), max_cutoff(), max_range(), and set_verlet_skin().

◆ update_ghosts_and_resort_particle()

void CellStructure::update_ghosts_and_resort_particle ( unsigned  data_parts)

Update ghost particles, with particle resort if needed.

Update ghost particles with data from the real particles. Resort particles if a resort is due.

Parameters
data_partsParticle parts to update, combination of Cells::DataPart

Definition at line 472 of file CellStructure.cpp.

References clear_resort_particles(), comm_cart, Cells::DATA_PART_BONDS, Cells::DATA_PART_PROPERTIES, get_local_particle(), ghost_particles(), ghosts_count(), ghosts_update(), Cells::RESORT_GLOBAL, Cells::RESORT_NONE, resort_particles(), and update_particle_index().

Referenced by correct_position_shake(), lb_tracers_add_particle_force_to_fluid(), and ScriptInterface::System::rotate_system().

◆ update_particle_index() [1/3]

void CellStructure::update_particle_index ( int  id,
Particle p 
)
inline

Update local particle index.

Update the entry for a particle in the local particle index.

Parameters
idEntry to update.
pPointer to the particle.

Definition at line 235 of file CellStructure.hpp.

References Particle::id().

Referenced by anonymous_namespace{CellStructure.cpp}::UpdateParticleIndexVisitor::operator()(), anonymous_namespace{CellStructure.cpp}::UpdateParticleIndexVisitor::operator()(), remove_particle(), update_ghosts_and_resort_particle(), update_particle_index(), and update_particle_index().

◆ update_particle_index() [2/3]

void CellStructure::update_particle_index ( Particle p)
inline

Update local particle index.

Update the entry for a particle in the local particle index.

Parameters
pPointer to the particle.

Definition at line 254 of file CellStructure.hpp.

References Particle::id(), and update_particle_index().

◆ update_particle_index() [3/3]

void CellStructure::update_particle_index ( ParticleList pl)
inline

Update local particle index.

Parameters
plList of particles whose index entries should be updated.

Definition at line 263 of file CellStructure.hpp.

References update_particle_index().

◆ update_verlet_stats()

void CellStructure::update_verlet_stats ( int  n_steps,
int  n_verlet_updates 
)
inline

Definition at line 577 of file CellStructure.hpp.

◆ use_parallel_for_each_local_particle()

bool CellStructure::use_parallel_for_each_local_particle ( ) const
inline

whether to use parallel version of for_each_local_particle

Definition at line 355 of file CellStructure.hpp.

Referenced by enumerate_local_particles(), and for_each_local_particle().

Member Data Documentation

◆ use_verlet_list

bool CellStructure::use_verlet_list = true

Definition at line 224 of file CellStructure.hpp.

Referenced by cabana_short_range(), non_bonded_loop(), and update_cabana_state().

◆ vector_length

constexpr auto CellStructure::vector_length = 1
staticconstexpr

Definition at line 172 of file CellStructure.hpp.


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