ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
Loading...
Searching...
No Matches
CellStructure.cpp
Go to the documentation of this file.
1/*
2 * Copyright (C) 2010-2026 The ESPResSo project
3 * Copyright (C) 2002,2003,2004,2005,2006,2007,2008,2009,2010
4 * Max-Planck-Institute for Polymer Research, Theory Group
5 *
6 * This file is part of ESPResSo.
7 *
8 * ESPResSo is free software: you can redistribute it and/or modify
9 * it under the terms of the GNU General Public License as published by
10 * the Free Software Foundation, either version 3 of the License, or
11 * (at your option) any later version.
12 *
13 * ESPResSo is distributed in the hope that it will be useful,
14 * but WITHOUT ANY WARRANTY; without even the implied warranty of
15 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16 * GNU General Public License for more details.
17 *
18 * You should have received a copy of the GNU General Public License
19 * along with this program. If not, see <http://www.gnu.org/licenses/>.
20 */
21
23
28
29#include "BoxGeometry.hpp"
30#include "LocalBox.hpp"
31#include "Particle.hpp"
32#include "aosoa_pack.hpp"
34#include "communication.hpp"
36#include "ghosts.hpp"
41#include "system/System.hpp"
42
43#include <utils/Vector.hpp>
44#include <utils/contains.hpp>
46#include <utils/math/sqr.hpp>
47
48#ifdef ESPRESSO_CALIPER
49#include <caliper/cali.h>
50#endif
51
52#include <boost/mpi/collectives/all_reduce.hpp>
53
54#ifdef ESPRESSO_SHARED_MEMORY_PARALLELISM
55#include <Cabana_Core.hpp>
56#include <Cabana_NeighborList.hpp>
57#include <Kokkos_Core.hpp>
58#include <omp.h>
59#endif
60
61#include <algorithm>
62#include <cassert>
63#include <cmath>
64#include <cstddef>
65#include <cstdint>
66#include <iterator>
67#include <memory>
68#include <numbers>
69#include <optional>
70#include <ranges>
71#include <set>
72#include <stdexcept>
73#include <string>
74#include <unordered_set>
75#include <utility>
76#include <variant>
77#include <vector>
78
80#ifdef ESPRESSO_SHARED_MEMORY_PARALLELISM
82 // Kokkos handle can only be freed after all Cabana containers have been freed
83 m_kokkos_handle.reset();
84#endif
85}
86
87#ifdef ESPRESSO_SHARED_MEMORY_PARALLELISM
89 m_local_force.reset();
90#ifdef ESPRESSO_ROTATION
91 m_local_torque.reset();
92#endif
93#ifdef ESPRESSO_NPT
94 m_local_virial.reset();
95#endif
96 m_id_to_index.reset();
97 m_aosoa.reset();
98 m_verlet_list_cabana.reset();
99 m_rebuild_verlet_list_cabana = true;
100}
101
102void CellStructure::set_kokkos_handle(std::shared_ptr<KokkosHandle> handle) {
103 m_kokkos_handle = std::move(handle);
104}
105
107 std::size_t number_of_unique_particles,
108 double local_box_volume,
109 std::size_t num_local_particles) {
110 if (std::isinf(pair_cutoff)) {
112 }
113 if (pair_cutoff < 0.) {
114 pair_cutoff = 0.;
115 }
116 // Estimate number of neighbors based on local density and cutoff sphere:
117 // volume n_neighbors = rho * (4/3) * pi * r^3, where rho = n_particles /
118 // volume
119 auto const local_density =
121 ? static_cast<double>(num_local_particles) / local_box_volume
122 : 0.;
123 auto const cutoff_sphere_volume =
124 (4. / 3.) * std::numbers::pi * Utils::int_pow<3>(pair_cutoff);
125 // account for local fluctuations. Empirical.
126 auto const fluctuation_factor = 2.;
127 auto max_counts = static_cast<std::size_t>(
129 std::size_t constexpr threshold_num = 16;
132 }
133 return max_counts;
134}
135
137#ifdef ESPRESSO_CALIPER
139#endif
140 assert(m_kokkos_handle);
141 using execution_space = Kokkos::DefaultExecutionSpace;
142 auto const num_threads = execution_space().concurrency();
143 auto const num_part = get_unique_particles().size();
144 auto const &system = get_system();
145 auto const local_box_volume = system.local_geo->volume();
148#ifdef ESPRESSO_COLLISION_DETECTION
149 if (system.has_collision_detection_enabled()) {
150 // TODO: use other types of Verlet list data structures
151 max_counts = num_part * 2ul;
152 }
153#endif
154 if (m_local_force) { // local properties are reallocated
155 Kokkos::realloc(get_local_force(), num_part, num_threads);
156#ifdef ESPRESSO_ROTATION
157 Kokkos::realloc(get_local_torque(), num_part, num_threads);
158#endif
159 Kokkos::realloc(get_id_to_index(), get_cached_max_local_particle_id() + 1);
160 Kokkos::deep_copy(get_id_to_index(), -1);
161 // Resize particle views using AoSoA_pack's resize method
162 m_aosoa->resize(num_part);
163 Kokkos::deep_copy(m_aosoa->flags, uint8_t{0});
164 m_verlet_list_cabana->reallocData(num_part, max_counts);
165 } else { // local properties are initialized
166 m_local_force =
167 std::make_unique<ForceType>("local_force", num_part, num_threads);
168#ifdef ESPRESSO_ROTATION
169 m_local_torque =
170 std::make_unique<ForceType>("local_torque", num_part, num_threads);
171#endif
172 m_id_to_index = std::make_unique<Kokkos::View<int *>>(
173 Kokkos::ViewAllocateWithoutInitializing("id_to_index"),
175 Kokkos::deep_copy(get_id_to_index(), -1);
176 // Create AoSoA_pack and initialize with resize
177 m_aosoa = std::make_unique<AoSoA_pack>();
178 m_aosoa->resize(num_part);
179 Kokkos::deep_copy(m_aosoa->flags, uint8_t{0});
180
181 m_verlet_list_cabana =
182 std::make_unique<ListType>(0ul, num_part, max_counts);
183 }
184#ifdef ESPRESSO_NPT
185 m_local_virial = std::make_unique<VirialType>("local_virial", num_threads);
186#endif
187}
188
190#ifdef ESPRESSO_CALIPER
192#endif
193 Kokkos::deep_copy(get_local_force(), 0.);
194}
195
197 Kokkos::deep_copy(get_local_force(), 0.);
198#ifdef ESPRESSO_ROTATION
199 Kokkos::deep_copy(get_local_torque(), 0.);
200#endif
201#ifdef ESPRESSO_NPT
202 Kokkos::deep_copy(get_local_virial(), 0.);
203#endif
204 Kokkos::deep_copy(get_aosoa().flags, uint8_t{0});
205}
206
208#ifdef ESPRESSO_CALIPER
210#endif
211 auto &unique_particles = m_unique_particles;
212 unique_particles.clear();
213 unique_particles.resize(count_local_particles());
214 std::unordered_set<int> registered_index{};
215 using execution_space = Kokkos::DefaultExecutionSpace;
216 int n_threads = execution_space().concurrency();
217 std::vector<int> max_ids(n_threads);
219 *this, [&unique_particles, &max_ids](std::size_t index, Particle &p) {
220 unique_particles[index] = &p;
221 const int thread_num = omp_get_thread_num();
222 max_ids[thread_num] = std::max(p.id(), max_ids[thread_num]);
223 });
224 int max_id = *(std::max_element(max_ids.begin(), max_ids.end()));
225 for (auto &p : ghost_particles()) {
226 auto const *local_particle = get_local_particle(p.id());
227 if (not local_particle) {
228 continue;
229 }
230 if (not local_particle->is_ghost()) {
231 continue;
232 }
233 if (registered_index.contains(p.id())) {
234 continue;
235 }
236 registered_index.insert(p.id());
237 unique_particles.emplace_back(&p);
238 max_id = std::max(p.id(), max_id);
239 }
240 registered_index.clear();
241 m_cached_max_local_particle_id = max_id;
242 m_num_local_particles_cached = unique_particles.size();
243}
244
245#endif // ESPRESSO_SHARED_MEMORY_PARALLELISM
246
248 : m_decomposition{std::make_unique<AtomDecomposition>(box)} {}
249
251 auto const max_id = get_max_local_particle_id();
252
253 for (auto const &p : local_particles()) {
254 auto const id = p.id();
255
257 throw std::runtime_error("Particle id out of bounds.");
258 }
259
260 if (get_local_particle(id) != &p) {
261 throw std::runtime_error("Invalid local particle index entry.");
262 }
263 }
264
265 /* checks: local particle id */
266 std::size_t local_part_cnt = 0u;
267 for (int n = 0; n < get_max_local_particle_id() + 1; n++) {
268 if (get_local_particle(n) != nullptr) {
270 if (get_local_particle(n)->id() != n) {
271 throw std::runtime_error("local_particles part has corrupted id.");
272 }
273 }
274 }
275
276 if (local_part_cnt != local_particles().size()) {
277 throw std::runtime_error(
278 std::to_string(local_particles().size()) + " parts in cells but " +
279 std::to_string(local_part_cnt) + " parts in local_particles");
280 }
281}
282
284 for (auto cell : decomposition().local_cells()) {
285 for (auto const &p : cell->particles()) {
286 if (particle_to_cell(p) != cell) {
287 throw std::runtime_error("misplaced particle with id " +
288 std::to_string(p.id()));
289 }
290 }
291 }
292}
293
295 auto remove_all_bonds_to = [id](BondList &bl) {
296 for (auto it = bl.begin(); it != bl.end();) {
297 if (Utils::contains(it->partner_ids(), id)) {
298 it = bl.erase(it);
299 } else {
300 std::advance(it, 1);
301 }
302 }
303 };
304
305 for (auto cell : decomposition().local_cells()) {
306 auto &parts = cell->particles();
307 for (auto it = parts.begin(); it != parts.end();) {
308 if (it->id() == id) {
309 it = parts.erase(it);
310 update_particle_index(id, nullptr);
312 } else {
313 remove_all_bonds_to(it->bonds());
314 it++;
315 }
316 }
317 }
318}
319
321 auto const sort_cell = particle_to_cell(p);
322 if (sort_cell) {
323 return std::addressof(
324 append_indexed_particle(sort_cell->particles(), std::move(p)));
325 }
326
327 return {};
328}
329
331 auto const sort_cell = particle_to_cell(p);
332 /* There is always at least one cell, so if the particle
333 * does not belong to a cell on this node we can put it there. */
334 auto cell = sort_cell ? sort_cell : decomposition().local_cells()[0];
335
336 /* If the particle isn't local a global resort may be
337 * needed, otherwise a local resort if sufficient. */
339
340 return std::addressof(
341 append_indexed_particle(cell->particles(), std::move(p)));
342}
343
345 auto it = std::ranges::find_if(std::ranges::views::reverse(m_particle_index),
346 [](auto const *p) { return p != nullptr; });
347
348 return (it != m_particle_index.rend()) ? (*it)->id() : -1;
349}
350
352 for (auto cell : decomposition().local_cells()) {
353 cell->particles().clear();
354 }
355
356 m_particle_index.clear();
357}
358
359/* Map the data parts flags from cells to those used internally
360 * by the ghost communication */
361unsigned map_data_parts(unsigned data_parts) {
362 using namespace Cells;
363
364 /* clang-format off */
365 return GHOSTTRANS_NONE
366 | ((data_parts & DATA_PART_PROPERTIES) ? GHOSTTRANS_PROPRTS : 0u)
367 | ((data_parts & DATA_PART_POSITION) ? GHOSTTRANS_POSITION : 0u)
368 | ((data_parts & DATA_PART_MOMENTUM) ? GHOSTTRANS_MOMENTUM : 0u)
369 | ((data_parts & DATA_PART_FORCE) ? GHOSTTRANS_FORCE : 0u)
371 | ((data_parts & DATA_PART_RATTLE) ? GHOSTTRANS_RATTLE : 0u)
372#endif
373 | ((data_parts & DATA_PART_BONDS) ? GHOSTTRANS_BONDS : 0u);
374 /* clang-format on */
375}
376
378 ghost_communicator(decomposition().exchange_ghosts_comm(),
379 *get_system().box_geo, GHOSTTRANS_PARTNUM);
380}
382 ghost_communicator(decomposition().exchange_ghosts_comm(),
384}
386 ghost_communicator(decomposition().collect_ghost_force_comm(),
387 *get_system().box_geo, GHOSTTRANS_FORCE);
388}
389#ifdef ESPRESSO_BOND_CONSTRAINT
394#endif
395
396namespace {
397/**
398 * @brief Apply a @ref ParticleChange to a particle index.
399 */
408} // namespace
409
411 invalidate_ghosts();
412
413 std::vector<ParticleChange> diff;
414
415 m_decomposition->resort(global_flag, diff);
416
417 for (auto d : diff) {
418 std::visit(UpdateParticleIndexVisitor{this}, d);
419 }
420
421 auto const &lebc = get_system().box_geo->lees_edwards_bc();
422 m_rebuild_verlet_list = true;
423 m_rebuild_verlet_list_cabana = true;
424 m_le_pos_offset_at_last_resort = lebc.pos_offset;
425
426#ifdef ESPRESSO_ADDITIONAL_CHECKS
429#endif
430}
431
433 auto &system = get_system();
434 auto &local_geo = *system.local_geo;
435 auto const &box_geo = *system.box_geo;
436 set_particle_decomposition(
437 std::make_unique<AtomDecomposition>(::comm_cart, box_geo));
439 local_geo.set_cell_structure_type(m_type);
440 system.on_cell_structure_change();
441}
442
444 double range, std::optional<std::pair<int, int>> fully_connected_boundary) {
445 auto &system = get_system();
446 auto &local_geo = *system.local_geo;
447 auto const &box_geo = *system.box_geo;
448 set_particle_decomposition(std::make_unique<RegularDecomposition>(
449 ::comm_cart, range, box_geo, local_geo, fully_connected_boundary));
451 local_geo.set_cell_structure_type(m_type);
452 system.on_cell_structure_change();
453}
454
456 std::set<int> n_square_types) {
457 auto &system = get_system();
458 auto &local_geo = *system.local_geo;
459 auto const &box_geo = *system.box_geo;
460 set_particle_decomposition(std::make_unique<HybridDecomposition>(
461 ::comm_cart, cutoff_regular, m_verlet_skin,
462 [&system]() { return system.get_global_ghost_flags(); }, box_geo,
463 local_geo, n_square_types));
465 local_geo.set_cell_structure_type(m_type);
466 system.on_cell_structure_change();
467}
468
470 assert(value >= 0.);
471 m_verlet_skin = value;
472 m_verlet_skin_set = true;
473 m_rebuild_verlet_list_cabana = true;
474 get_system().on_verlet_skin_change();
475}
476
479 auto const max_cut = get_system().maximal_cutoff();
480 if (max_cut <= 0.) {
481 throw std::runtime_error(
482 "cannot automatically determine skin, please set it manually");
483 }
484 /* maximal skin that can be used without resorting is the maximal
485 * range of the cell system minus what is needed for interactions. */
486 auto const max_range = std::ranges::min(max_cutoff());
487 auto const new_skin = std::min(0.4 * max_cut, max_range - max_cut);
489}
490
492 /* data parts that are only updated on resort */
493 auto constexpr resort_only_parts =
495
496 auto const global_resort = boost::mpi::all_reduce(
497 ::comm_cart, m_resort_particles, std::bit_or<unsigned>());
498
501
502 /* Resort cell system */
504 ghosts_count();
506
507 /* Add the ghost particles to the index if we don't already
508 * have them. */
509 for (auto &p : ghost_particles()) {
510 if (get_local_particle(p.id()) == nullptr) {
511 update_particle_index(p.id(), &p);
512 }
513 }
514
515 /* Particles are now sorted */
517 } else {
518 /* Communication step: ghost information */
520 }
521}
522
523#ifdef ESPRESSO_SHARED_MEMORY_PARALLELISM
524void CellStructure::parallel_for_each_particle_impl(
525 std::span<Cell *const> cells, ParticleUnaryOp &f) const {
526 if (cells.size() > 1) {
527 Kokkos::parallel_for( // loop over cells
528 "for_each_local_particle", cells.size(), [&](auto cell_idx) {
529 for (auto &p : cells[cell_idx]->particles())
530 f(p);
531 });
532 } else if (cells.size() == 1) {
533 auto &particles = cells.front()->particles();
534 Kokkos::parallel_for( // loop over particles
535 "for_each_local_particle", particles.size(),
536 [&](auto part_idx) { f(*(particles.begin() + part_idx)); });
537 }
538}
539#endif // ESPRESSO_SHARED_MEMORY_PARALLELISM
540
542 Utils::Vector3d const &additional_offset) const {
543 auto const lim = Utils::sqr(m_verlet_skin / 2.) - additional_offset.norm2();
544
546 [lim](bool &result, Particle const &p) {
547 if ((p.pos() - p.pos_at_last_verlet_update()).norm2() > lim) {
548 result = true;
549 }
550 };
551
552 Reduction::ReductionOp<bool> reduce_op = [](bool &acc, bool const &val) {
553 acc |= val;
554 };
555
557}
@ 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.
Bond storage.
Definition BondList.hpp:84
Describes a cell structure / cell system.
ParticleRange ghost_particles() const
auto & get_local_force()
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.
auto & get_id_to_index()
std::size_t count_local_particles() const
virtual ~CellStructure()
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
cudaStream_t stream[1]
CUDA streams for parallel computing on CPU and GPU.
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.
Definition ghosts.cpp:443
Ghost particles and particle exchange.
@ GHOSTTRANS_MOMENTUM
transfer ParticleMomentum
Definition ghosts.hpp:132
@ GHOSTTRANS_RATTLE
transfer ParticleRattle
Definition ghosts.hpp:137
@ GHOSTTRANS_PARTNUM
resize the receiver particle arrays to the size of the senders
Definition ghosts.hpp:140
@ GHOSTTRANS_POSITION
transfer ParticlePosition
Definition ghosts.hpp:130
@ GHOSTTRANS_PROPRTS
transfer ParticleProperties
Definition ghosts.hpp:128
@ GHOSTTRANS_FORCE
transfer ParticleForce
Definition ghosts.hpp:134
@ GHOSTTRANS_NONE
Definition ghosts.hpp:126
@ GHOSTTRANS_BONDS
Definition ghosts.hpp:141
@ 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.
Definition sqr.hpp:28
bool contains(Range &&rng, T const &value)
Check whether a range contains a value.
Definition contains.hpp:36
STL namespace.
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.
Definition Particle.hpp:435
auto const & id() const
Definition Particle.hpp:454