34#include "communication.hpp"
40#include "system/System.hpp"
47#ifdef ESPRESSO_CALIPER
48#include <caliper/cali.h>
51#include <boost/mpi/collectives/all_reduce.hpp>
53#ifdef ESPRESSO_SHARED_MEMORY_PARALLELISM
54#include <Cabana_Core.hpp>
55#include <Cabana_NeighborList.hpp>
56#include <Kokkos_Core.hpp>
78#ifdef ESPRESSO_SHARED_MEMORY_PARALLELISM
81 m_kokkos_handle.reset();
85#ifdef ESPRESSO_SHARED_MEMORY_PARALLELISM
87 m_local_force.reset();
88#ifdef ESPRESSO_ROTATION
89 m_local_torque.reset();
92 m_local_virial.reset();
94 m_id_to_index.reset();
96 m_verlet_list_cabana.reset();
97 m_rebuild_verlet_list_cabana =
true;
101 m_kokkos_handle = std::move(handle);
105 std::size_t number_of_unique_particles,
106 double local_box_volume,
107 std::size_t num_local_particles) {
108 if (std::isinf(pair_cutoff)) {
109 return number_of_unique_particles;
111 if (pair_cutoff < 0.) {
117 auto const local_density =
118 (local_box_volume > 0. && num_local_particles > 0)
119 ?
static_cast<double>(num_local_particles) / local_box_volume
121 auto const cutoff_sphere_volume =
122 (4. / 3.) * std::numbers::pi * Utils::int_pow<3>(pair_cutoff);
124 auto const fluctuation_factor = 2.;
125 auto max_counts =
static_cast<std::size_t
>(
126 std::ceil(fluctuation_factor * local_density * cutoff_sphere_volume));
127 std::size_t
constexpr threshold_num = 16;
128 if (max_counts < threshold_num) {
129 max_counts = std::min(threshold_num, number_of_unique_particles);
135#ifdef ESPRESSO_CALIPER
136 CALI_CXX_MARK_FUNCTION;
138 assert(m_kokkos_handle);
139 using execution_space = Kokkos::DefaultExecutionSpace;
140 auto const num_threads = execution_space().concurrency();
143 auto const local_box_volume = system.local_geo->volume();
146#ifdef ESPRESSO_COLLISION_DETECTION
147 if (system.has_collision_detection_enabled()) {
149 max_counts = num_part * 2ul;
154#ifdef ESPRESSO_ROTATION
160 m_aosoa->resize(num_part);
161 Kokkos::deep_copy(m_aosoa->flags, uint8_t{0});
162 m_verlet_list_cabana->reallocData(num_part, max_counts);
165 std::make_unique<ForceType>(
"local_force", num_part, num_threads);
166#ifdef ESPRESSO_ROTATION
168 std::make_unique<ForceType>(
"local_torque", num_part, num_threads);
170 m_id_to_index = std::make_unique<Kokkos::View<int *>>(
171 Kokkos::ViewAllocateWithoutInitializing(
"id_to_index"),
175 m_aosoa = std::make_unique<AoSoA_pack>();
176 m_aosoa->resize(num_part);
177 Kokkos::deep_copy(m_aosoa->flags, uint8_t{0});
179 m_verlet_list_cabana =
180 std::make_unique<ListType>(0ul, num_part, max_counts);
183 m_local_virial = std::make_unique<VirialType>(
"local_virial", num_threads);
188#ifdef ESPRESSO_CALIPER
189 CALI_CXX_MARK_FUNCTION;
196#ifdef ESPRESSO_ROTATION
202 Kokkos::deep_copy(
get_aosoa().flags, uint8_t{0});
206#ifdef ESPRESSO_CALIPER
207 CALI_CXX_MARK_FUNCTION;
209 auto &unique_particles = m_unique_particles;
210 unique_particles.clear();
212 std::unordered_set<int> registered_index{};
213 using execution_space = Kokkos::DefaultExecutionSpace;
214 int n_threads = execution_space().concurrency();
215 std::vector<int> max_ids(n_threads);
217 *
this, [&unique_particles, &max_ids](std::size_t index,
Particle &p) {
218 unique_particles[index] = &p;
219 const int thread_num = omp_get_thread_num();
220 max_ids[thread_num] = std::max(p.
id(), max_ids[thread_num]);
222 int max_id = *(std::max_element(max_ids.begin(), max_ids.end()));
225 if (not local_particle) {
228 if (not local_particle->is_ghost()) {
231 if (registered_index.contains(p.
id())) {
234 registered_index.insert(p.
id());
235 unique_particles.emplace_back(&p);
236 max_id = std::max(p.
id(), max_id);
238 registered_index.clear();
239 m_cached_max_local_particle_id = max_id;
240 m_num_local_particles_cached = unique_particles.size();
252 auto const id = p.id();
254 if (id < 0 || id > max_id) {
255 throw std::runtime_error(
"Particle id out of bounds.");
259 throw std::runtime_error(
"Invalid local particle index entry.");
264 std::size_t local_part_cnt = 0u;
269 throw std::runtime_error(
"local_particles part has corrupted id.");
275 throw std::runtime_error(
277 std::to_string(local_part_cnt) +
" parts in local_particles");
283 for (
auto const &p : cell->particles()) {
284 if (particle_to_cell(p) != cell) {
285 throw std::runtime_error(
"misplaced particle with id " +
286 std::to_string(p.id()));
293 auto remove_all_bonds_to = [id](
BondList &bl) {
294 for (
auto it = bl.begin(); it != bl.end();) {
304 auto &parts = cell->particles();
305 for (
auto it = parts.begin(); it != parts.end();) {
306 if (it->id() == id) {
307 it = parts.erase(it);
311 remove_all_bonds_to(it->bonds());
319 auto const sort_cell = particle_to_cell(p);
321 return std::addressof(
322 append_indexed_particle(sort_cell->particles(), std::move(p)));
329 auto const sort_cell = particle_to_cell(p);
338 return std::addressof(
339 append_indexed_particle(cell->particles(), std::move(p)));
343 auto it = std::ranges::find_if(std::ranges::views::reverse(m_particle_index),
344 [](
auto const *p) {
return p !=
nullptr; });
346 return (it != m_particle_index.rend()) ? (*it)->id() : -1;
351 cell->particles().clear();
354 m_particle_index.clear();
360 using namespace Cells;
368#ifdef ESPRESSO_BOND_CONSTRAINT
387#ifdef ESPRESSO_BOND_CONSTRAINT
411 std::vector<ParticleChange> diff;
413 m_decomposition->resort(global_flag, diff);
415 for (
auto d : diff) {
416 std::visit(UpdateParticleIndexVisitor{
this}, d);
419 auto const &lebc =
get_system().box_geo->lees_edwards_bc();
420 m_rebuild_verlet_list =
true;
421 m_rebuild_verlet_list_cabana =
true;
422 m_le_pos_offset_at_last_resort = lebc.pos_offset;
424#ifdef ESPRESSO_ADDITIONAL_CHECKS
432 auto &local_geo = *system.local_geo;
433 auto const &box_geo = *system.box_geo;
434 set_particle_decomposition(
435 std::make_unique<AtomDecomposition>(
::comm_cart, box_geo));
437 local_geo.set_cell_structure_type(m_type);
438 system.on_cell_structure_change();
442 double range, std::optional<std::pair<int, int>> fully_connected_boundary) {
444 auto &local_geo = *system.local_geo;
445 auto const &box_geo = *system.box_geo;
446 set_particle_decomposition(std::make_unique<RegularDecomposition>(
447 ::comm_cart, range, box_geo, local_geo, fully_connected_boundary));
449 local_geo.set_cell_structure_type(m_type);
450 system.on_cell_structure_change();
454 std::set<int> n_square_types) {
456 auto &local_geo = *system.local_geo;
457 auto const &box_geo = *system.box_geo;
458 set_particle_decomposition(std::make_unique<HybridDecomposition>(
460 [&system]() {
return system.get_global_ghost_flags(); }, box_geo,
461 local_geo, n_square_types));
463 local_geo.set_cell_structure_type(m_type);
464 system.on_cell_structure_change();
469 m_verlet_skin = value;
470 m_verlet_skin_set =
true;
471 m_rebuild_verlet_list_cabana =
true;
477 auto const max_cut =
get_system().maximal_cutoff();
479 throw std::runtime_error(
480 "cannot automatically determine skin, please set it manually");
485 auto const new_skin = std::min(0.4 * max_cut,
max_range - max_cut);
491 auto constexpr resort_only_parts =
494 auto const global_resort = boost::mpi::all_reduce(
495 ::comm_cart, m_resort_particles, std::bit_or<unsigned>());
521#ifdef ESPRESSO_SHARED_MEMORY_PARALLELISM
522void CellStructure::parallel_for_each_particle_impl(
524 if (cells.size() > 1) {
525 Kokkos::parallel_for(
526 "for_each_local_particle", cells.size(), [&](
auto cell_idx) {
527 for (auto &p : cells[cell_idx]->particles())
530 }
else if (cells.size() == 1) {
531 auto &
particles = cells.front()->particles();
532 Kokkos::parallel_for(
534 [&](
auto part_idx) { f(*(particles.begin() + part_idx)); });
541 auto const lim =
Utils::sqr(m_verlet_skin / 2.) - additional_offset.
norm2();
544 [lim](
bool &result,
Particle const &p) {
545 if ((p.pos() - p.pos_at_last_verlet_update()).norm2() > lim) {
@ NSQUARE
Atom decomposition (N-square).
@ HYBRID
Hybrid decomposition.
@ REGULAR
Regular decomposition.
unsigned map_data_parts(unsigned data_parts)
Map the data parts flags from cells to those used internally by the ghost communication.
static auto estimate_max_counts(double pair_cutoff, std::size_t number_of_unique_particles, double local_box_volume, std::size_t num_local_particles)
unsigned map_data_parts(unsigned data_parts)
Map the data parts flags from cells to those used internally by the ghost communication.
std::function< void(Particle &)> ParticleUnaryOp
Vector implementation and trait types for boost qvm interoperability.
Atom decomposition cell system.
Describes a cell structure / cell system.
ParticleRange ghost_particles() const
Particle * get_local_particle(int id)
Get a local particle by id.
void set_kokkos_handle(std::shared_ptr< KokkosHandle > handle)
void check_particle_sorting() const
Check that particles are in the correct cell.
std::size_t count_local_particles() const
void clear_resort_particles()
Set the resort level to sorted.
auto is_verlet_skin_set() const
Whether the Verlet skin is set.
void clear_local_properties()
ParticleDecomposition const & decomposition() const
Get the underlying particle decomposition.
void update_ghosts_and_resort_particle(unsigned data_parts)
Update ghost particles, with particle resort if needed.
Particle * add_local_particle(Particle &&p)
Add a particle.
void set_verlet_skin_heuristic()
Set the Verlet skin using a heuristic.
void set_verlet_skin(double value)
Set the Verlet skin.
void ghosts_update(unsigned data_parts)
Update ghost particles.
int get_cached_max_local_particle_id() const
CellStructure(BoxGeometry const &box)
auto & get_local_torque()
auto & get_local_virial()
void update_particle_index(int id, Particle *p)
Update local particle index.
void ghosts_reduce_forces()
Add forces and torques from ghost particles to real particles.
auto const & get_unique_particles() const
void rebuild_local_properties(double pair_cutoff)
Utils::Vector3d max_range() const
Maximal pair range supported by current cell system.
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,...
void ghosts_count()
Synchronize number of ghosts.
void set_resort_particles(Cells::Resort level)
Increase the local resort level at least to level.
void remove_particle(int id)
Remove a particle.
Particle * add_particle(Particle &&p)
Add a particle.
std::size_t get_num_local_particles_cached() const
void resort_particles(bool global_flag)
Resort particles.
void check_particle_index() const
Check that particle index is commensurate with particles.
void set_regular_decomposition(double range, std::optional< std::pair< int, int > > fully_connected_boundary)
Set the particle decomposition to RegularDecomposition.
void set_atom_decomposition()
Set the particle decomposition to AtomDecomposition.
void remove_all_particles()
Remove all particles from the cell system.
ParticleRange local_particles() const
void ghosts_reduce_rattle_correction()
Add rattle corrections from ghost particles to real particles.
void set_hybrid_decomposition(double cutoff_regular, std::set< int > n_square_types)
Set the particle decomposition to HybridDecomposition.
int get_max_local_particle_id() const
Get the maximal particle id on this node.
Utils::Vector3d max_cutoff() const
Maximal cutoff supported by current cell system.
void reset_local_properties()
virtual std::span< Cell *const > local_cells() const =0
Get pointer to local cells.
base_type::size_type size() const
constexpr T norm2() const
boost::mpi::communicator comm_cart
The communicator.
void ghost_communicator(GhostCommunicator const &gcr, BoxGeometry const &box_geo, unsigned int data_parts)
Do a ghost communication with the specified data parts.
@ GHOSTTRANS_MOMENTUM
transfer ParticleMomentum
@ GHOSTTRANS_RATTLE
transfer ParticleRattle
@ GHOSTTRANS_PARTNUM
resize the receiver particle arrays to the size of the senders
@ GHOSTTRANS_POSITION
transfer ParticlePosition
@ GHOSTTRANS_PROPRTS
transfer ParticleProperties
@ GHOSTTRANS_FORCE
transfer ParticleForce
@ DATA_PART_PROPERTIES
Particle::p.
@ DATA_PART_BONDS
Particle::bonds.
ParticleRange particles(std::span< Cell *const > cells)
std::function< void(ResultType &, ResultType const &)> ReductionOp
Join two partial reduction results.
std::function< void(ResultType &, Particle const &)> AddPartialResultKernel
Kernel that adds the result from a single particle to a reduction.
DEVICE_QUALIFIER constexpr T sqr(T x)
Calculates the SQuaRe of x.
bool contains(Range &&rng, T const &value)
Check whether a range contains a value.
void enumerate_local_particles(CellStructure const &cs, Kernel &&kernel)
Run a kernel on all local particles with enumeration.
ResultType reduce_over_local_particles(CellStructure const &cs, Reduction::AddPartialResultKernel< ResultType > add_partial, Reduction::ReductionOp< ResultType > reduce_op)
performs a reduction over all particles
Struct holding all information for one particle.
Apply a ParticleChange to a particle index.
void operator()(RemovedParticle rp) const
void operator()(ModifiedList mp) const