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-2022 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"
40#include "system/System.hpp"
41
42#include <utils/Vector.hpp>
43#include <utils/contains.hpp>
45#include <utils/math/sqr.hpp>
46
47#ifdef ESPRESSO_CALIPER
48#include <caliper/cali.h>
49#endif
50
51#include <boost/mpi/collectives/all_reduce.hpp>
52
53#ifdef ESPRESSO_SHARED_MEMORY_PARALLELISM
54#include <Cabana_Core.hpp>
55#include <Cabana_NeighborList.hpp>
56#include <Kokkos_Core.hpp>
57#include <omp.h>
58#endif
59
60#include <algorithm>
61#include <cassert>
62#include <cmath>
63#include <cstddef>
64#include <cstdint>
65#include <iterator>
66#include <memory>
67#include <numbers>
68#include <optional>
69#include <ranges>
70#include <set>
71#include <stdexcept>
72#include <string>
73#include <utility>
74#include <variant>
75#include <vector>
76
78#ifdef ESPRESSO_SHARED_MEMORY_PARALLELISM
80 // Kokkos handle can only be freed after all Cabana containers have been freed
81 m_kokkos_handle.reset();
82#endif
83}
84
85#ifdef ESPRESSO_SHARED_MEMORY_PARALLELISM
87 m_local_force.reset();
88#ifdef ESPRESSO_ROTATION
89 m_local_torque.reset();
90#endif
91#ifdef ESPRESSO_NPT
92 m_local_virial.reset();
93#endif
94 m_id_to_index.reset();
95 m_aosoa.reset();
96 m_verlet_list_cabana.reset();
97 m_rebuild_verlet_list_cabana = true;
98}
99
100void CellStructure::set_kokkos_handle(std::shared_ptr<KokkosHandle> handle) {
101 m_kokkos_handle = std::move(handle);
102}
103
104static auto estimate_max_counts(double pair_cutoff,
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;
110 }
111 if (pair_cutoff < 0.) {
112 pair_cutoff = 0.;
113 }
114 // Estimate number of neighbors based on local density and cutoff sphere:
115 // volume n_neighbors = rho * (4/3) * pi * r^3, where rho = n_particles /
116 // volume
117 auto const local_density =
118 (local_box_volume > 0. && num_local_particles > 0)
119 ? static_cast<double>(num_local_particles) / local_box_volume
120 : 0.;
121 auto const cutoff_sphere_volume =
122 (4. / 3.) * std::numbers::pi * Utils::int_pow<3>(pair_cutoff);
123 // account for local fluctuations. Empirical.
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);
130 }
131 return max_counts;
132}
133
134void CellStructure::rebuild_local_properties(double const pair_cutoff) {
135#ifdef ESPRESSO_CALIPER
136 CALI_CXX_MARK_FUNCTION;
137#endif
138 assert(m_kokkos_handle);
139 using execution_space = Kokkos::DefaultExecutionSpace;
140 auto const num_threads = execution_space().concurrency();
141 auto const num_part = get_unique_particles().size();
142 auto const &system = get_system();
143 auto const local_box_volume = system.local_geo->volume();
144 auto max_counts = estimate_max_counts(pair_cutoff, num_part, local_box_volume,
146#ifdef ESPRESSO_COLLISION_DETECTION
147 if (system.has_collision_detection_enabled()) {
148 // TODO: use other types of Verlet list data structures
149 max_counts = num_part * 2ul;
150 }
151#endif
152 if (m_local_force) { // local properties are reallocated
153 Kokkos::realloc(get_local_force(), num_part, num_threads);
154#ifdef ESPRESSO_ROTATION
155 Kokkos::realloc(get_local_torque(), num_part, num_threads);
156#endif
157 Kokkos::realloc(get_id_to_index(), get_cached_max_local_particle_id() + 1);
158 Kokkos::deep_copy(get_id_to_index(), -1);
159 // Resize particle views using AoSoA_pack's resize method
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);
163 } else { // local properties are initialized
164 m_local_force =
165 std::make_unique<ForceType>("local_force", num_part, num_threads);
166#ifdef ESPRESSO_ROTATION
167 m_local_torque =
168 std::make_unique<ForceType>("local_torque", num_part, num_threads);
169#endif
170 m_id_to_index = std::make_unique<Kokkos::View<int *>>(
171 Kokkos::ViewAllocateWithoutInitializing("id_to_index"),
173 Kokkos::deep_copy(get_id_to_index(), -1);
174 // Create AoSoA_pack and initialize with resize
175 m_aosoa = std::make_unique<AoSoA_pack>();
176 m_aosoa->resize(num_part);
177 Kokkos::deep_copy(m_aosoa->flags, uint8_t{0});
178
179 m_verlet_list_cabana =
180 std::make_unique<ListType>(0ul, num_part, max_counts);
181 }
182#ifdef ESPRESSO_NPT
183 m_local_virial = std::make_unique<VirialType>("local_virial", num_threads);
184#endif
185}
186
188#ifdef ESPRESSO_CALIPER
189 CALI_CXX_MARK_FUNCTION;
190#endif
191 Kokkos::deep_copy(get_local_force(), 0.);
192}
193
195 Kokkos::deep_copy(get_local_force(), 0.);
196#ifdef ESPRESSO_ROTATION
197 Kokkos::deep_copy(get_local_torque(), 0.);
198#endif
199#ifdef ESPRESSO_NPT
200 Kokkos::deep_copy(get_local_virial(), 0.);
201#endif
202 Kokkos::deep_copy(get_aosoa().flags, uint8_t{0});
203}
204
206#ifdef ESPRESSO_CALIPER
207 CALI_CXX_MARK_FUNCTION;
208#endif
209 auto &unique_particles = m_unique_particles;
210 unique_particles.clear();
211 unique_particles.resize(count_local_particles());
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]);
221 });
222 int max_id = *(std::max_element(max_ids.begin(), max_ids.end()));
223 for (auto &p : ghost_particles()) {
224 auto const *local_particle = get_local_particle(p.id());
225 if (not local_particle) {
226 continue;
227 }
228 if (not local_particle->is_ghost()) {
229 continue;
230 }
231 if (registered_index.contains(p.id())) {
232 continue;
233 }
234 registered_index.insert(p.id());
235 unique_particles.emplace_back(&p);
236 max_id = std::max(p.id(), max_id);
237 }
238 registered_index.clear();
239 m_cached_max_local_particle_id = max_id;
240 m_num_local_particles_cached = unique_particles.size();
241}
242
243#endif // ESPRESSO_SHARED_MEMORY_PARALLELISM
244
246 : m_decomposition{std::make_unique<AtomDecomposition>(box)} {}
247
249 auto const max_id = get_max_local_particle_id();
250
251 for (auto const &p : local_particles()) {
252 auto const id = p.id();
253
254 if (id < 0 || id > max_id) {
255 throw std::runtime_error("Particle id out of bounds.");
256 }
257
258 if (get_local_particle(id) != &p) {
259 throw std::runtime_error("Invalid local particle index entry.");
260 }
261 }
262
263 /* checks: local particle id */
264 std::size_t local_part_cnt = 0u;
265 for (int n = 0; n < get_max_local_particle_id() + 1; n++) {
266 if (get_local_particle(n) != nullptr) {
267 local_part_cnt++;
268 if (get_local_particle(n)->id() != n) {
269 throw std::runtime_error("local_particles part has corrupted id.");
270 }
271 }
272 }
273
274 if (local_part_cnt != local_particles().size()) {
275 throw std::runtime_error(
276 std::to_string(local_particles().size()) + " parts in cells but " +
277 std::to_string(local_part_cnt) + " parts in local_particles");
278 }
279}
280
282 for (auto cell : decomposition().local_cells()) {
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()));
287 }
288 }
289 }
290}
291
293 auto remove_all_bonds_to = [id](BondList &bl) {
294 for (auto it = bl.begin(); it != bl.end();) {
295 if (Utils::contains(it->partner_ids(), id)) {
296 it = bl.erase(it);
297 } else {
298 std::advance(it, 1);
299 }
300 }
301 };
302
303 for (auto cell : decomposition().local_cells()) {
304 auto &parts = cell->particles();
305 for (auto it = parts.begin(); it != parts.end();) {
306 if (it->id() == id) {
307 it = parts.erase(it);
308 update_particle_index(id, nullptr);
310 } else {
311 remove_all_bonds_to(it->bonds());
312 it++;
313 }
314 }
315 }
316}
317
319 auto const sort_cell = particle_to_cell(p);
320 if (sort_cell) {
321 return std::addressof(
322 append_indexed_particle(sort_cell->particles(), std::move(p)));
323 }
324
325 return {};
326}
327
329 auto const sort_cell = particle_to_cell(p);
330 /* There is always at least one cell, so if the particle
331 * does not belong to a cell on this node we can put it there. */
332 auto cell = sort_cell ? sort_cell : decomposition().local_cells()[0];
333
334 /* If the particle isn't local a global resort may be
335 * needed, otherwise a local resort if sufficient. */
337
338 return std::addressof(
339 append_indexed_particle(cell->particles(), std::move(p)));
340}
341
343 auto it = std::ranges::find_if(std::ranges::views::reverse(m_particle_index),
344 [](auto const *p) { return p != nullptr; });
345
346 return (it != m_particle_index.rend()) ? (*it)->id() : -1;
347}
348
350 for (auto cell : decomposition().local_cells()) {
351 cell->particles().clear();
352 }
353
354 m_particle_index.clear();
355}
356
357/* Map the data parts flags from cells to those used internally
358 * by the ghost communication */
359unsigned map_data_parts(unsigned data_parts) {
360 using namespace Cells;
361
362 /* clang-format off */
363 return GHOSTTRANS_NONE
364 | ((data_parts & DATA_PART_PROPERTIES) ? GHOSTTRANS_PROPRTS : 0u)
365 | ((data_parts & DATA_PART_POSITION) ? GHOSTTRANS_POSITION : 0u)
366 | ((data_parts & DATA_PART_MOMENTUM) ? GHOSTTRANS_MOMENTUM : 0u)
367 | ((data_parts & DATA_PART_FORCE) ? GHOSTTRANS_FORCE : 0u)
368#ifdef ESPRESSO_BOND_CONSTRAINT
369 | ((data_parts & DATA_PART_RATTLE) ? GHOSTTRANS_RATTLE : 0u)
370#endif
371 | ((data_parts & DATA_PART_BONDS) ? GHOSTTRANS_BONDS : 0u);
372 /* clang-format on */
373}
374
376 ghost_communicator(decomposition().exchange_ghosts_comm(),
377 *get_system().box_geo, GHOSTTRANS_PARTNUM);
378}
379void CellStructure::ghosts_update(unsigned data_parts) {
380 ghost_communicator(decomposition().exchange_ghosts_comm(),
381 *get_system().box_geo, map_data_parts(data_parts));
382}
384 ghost_communicator(decomposition().collect_ghost_force_comm(),
385 *get_system().box_geo, GHOSTTRANS_FORCE);
386}
387#ifdef ESPRESSO_BOND_CONSTRAINT
392#endif
393
394namespace {
395/**
396 * @brief Apply a @ref ParticleChange to a particle index.
397 */
400
402 cs->update_particle_index(rp.id, nullptr);
403 }
405};
406} // namespace
407
408void CellStructure::resort_particles(bool global_flag) {
409 invalidate_ghosts();
410
411 std::vector<ParticleChange> diff;
412
413 m_decomposition->resort(global_flag, diff);
414
415 for (auto d : diff) {
416 std::visit(UpdateParticleIndexVisitor{this}, d);
417 }
418
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;
423
424#ifdef ESPRESSO_ADDITIONAL_CHECKS
427#endif
428}
429
431 auto &system = get_system();
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();
439}
440
442 double range, std::optional<std::pair<int, int>> fully_connected_boundary) {
443 auto &system = get_system();
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();
451}
452
454 std::set<int> n_square_types) {
455 auto &system = get_system();
456 auto &local_geo = *system.local_geo;
457 auto const &box_geo = *system.box_geo;
458 set_particle_decomposition(std::make_unique<HybridDecomposition>(
459 ::comm_cart, cutoff_regular, m_verlet_skin,
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();
465}
466
468 assert(value >= 0.);
469 m_verlet_skin = value;
470 m_verlet_skin_set = true;
471 m_rebuild_verlet_list_cabana = true;
472 get_system().on_verlet_skin_change();
473}
474
476 assert(not is_verlet_skin_set());
477 auto const max_cut = get_system().maximal_cutoff();
478 if (max_cut <= 0.) {
479 throw std::runtime_error(
480 "cannot automatically determine skin, please set it manually");
481 }
482 /* maximal skin that can be used without resorting is the maximal
483 * range of the cell system minus what is needed for interactions. */
484 auto const max_range = std::ranges::min(max_cutoff());
485 auto const new_skin = std::min(0.4 * max_cut, max_range - max_cut);
486 set_verlet_skin(new_skin);
487}
488
490 /* data parts that are only updated on resort */
491 auto constexpr resort_only_parts =
493
494 auto const global_resort = boost::mpi::all_reduce(
495 ::comm_cart, m_resort_particles, std::bit_or<unsigned>());
496
497 if (global_resort != Cells::RESORT_NONE) {
498 auto const do_global_resort = (global_resort & Cells::RESORT_GLOBAL) != 0;
499
500 /* Resort cell system */
501 resort_particles(do_global_resort);
502 ghosts_count();
503 ghosts_update(data_parts);
504
505 /* Add the ghost particles to the index if we don't already
506 * have them. */
507 for (auto &p : ghost_particles()) {
508 if (get_local_particle(p.id()) == nullptr) {
509 update_particle_index(p.id(), &p);
510 }
511 }
512
513 /* Particles are now sorted */
515 } else {
516 /* Communication step: ghost information */
517 ghosts_update(data_parts & ~resort_only_parts);
518 }
519}
520
521#ifdef ESPRESSO_SHARED_MEMORY_PARALLELISM
522void CellStructure::parallel_for_each_particle_impl(
523 std::span<Cell *const> cells, ParticleUnaryOp &f) const {
524 if (cells.size() > 1) {
525 Kokkos::parallel_for( // loop over cells
526 "for_each_local_particle", cells.size(), [&](auto cell_idx) {
527 for (auto &p : cells[cell_idx]->particles())
528 f(p);
529 });
530 } else if (cells.size() == 1) {
531 auto &particles = cells.front()->particles();
532 Kokkos::parallel_for( // loop over particles
533 "for_each_local_particle", particles.size(),
534 [&](auto part_idx) { f(*(particles.begin() + part_idx)); });
535 }
536}
537#endif // ESPRESSO_SHARED_MEMORY_PARALLELISM
538
540 Utils::Vector3d const &additional_offset) const {
541 auto const lim = Utils::sqr(m_verlet_skin / 2.) - additional_offset.norm2();
542
544 [lim](bool &result, Particle const &p) {
545 if ((p.pos() - p.pos_at_last_verlet_update()).norm2() > lim) {
546 result = true;
547 }
548 };
549
550 Reduction::ReductionOp<bool> reduce_op = [](bool &acc, bool const &val) {
551 acc |= val;
552 };
553
554 return reduce_over_local_particles(*this, add_partial, reduce_op);
555}
@ 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
constexpr T norm2() const
Definition Vector.hpp:159
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:442
@ 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