ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
Loading...
Searching...
No Matches
CellStructure.hpp
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
22#pragma once
23
25
26#include "BoxGeometry.hpp"
27#include "LocalBox.hpp"
28#include "Particle.hpp"
29#include "ParticleList.hpp"
30#include "ParticleRange.hpp"
32#include "bond_error.hpp"
33#include "cell_system/Cell.hpp"
35#include "config/config.hpp"
36#include "ghosts.hpp"
37#include "system/Leaf.hpp"
38
39#include <utils/Vector.hpp>
40
41#include <boost/container/static_vector.hpp>
42#include <boost/iterator/indirect_iterator.hpp>
43#include <boost/range/algorithm/transform.hpp>
44
45#include <algorithm>
46#include <cassert>
47#include <concepts>
48#include <cstddef>
49#include <functional>
50#include <iterator>
51#include <memory>
52#include <optional>
53#include <set>
54#include <span>
55#include <stdexcept>
56#include <unordered_set>
57#include <utility>
58#include <vector>
59
60#ifdef ESPRESSO_CALIPER
61#include <caliper/cali.h>
62#endif
63
64// forward declarations
65namespace Kokkos {
66template <class DataType, class... Properties> class View;
67class HostSpace;
68struct LayoutRight;
69template <unsigned T> struct MemoryTraits;
70} // namespace Kokkos
71namespace Cabana {
72class HalfNeighborTag;
73struct VerletLayout2D;
74class TeamVectorOpTag;
75template <typename... Types> struct MemberTypes;
76template <class DataType, class MemorySpace, int, class MemoryTraits>
77class AoSoA;
78} // namespace Cabana
79struct KokkosHandle;
80template <class MemorySpace, class ListAlgorithm, class Layout, class BuildTag>
82struct LocalBondState;
83
84template <typename Callable>
85concept ParticleCallback = requires(Callable c, Particle &p) {
86 { c(p) } -> std::same_as<void>;
87};
88
89using ParticleUnaryOp = std::function<void(Particle &)>;
90
91namespace Cells {
92enum Resort : unsigned {
95 RESORT_GLOBAL = 2u
96};
97
98/**
99 * @brief Flags to select particle parts for communication.
100 */
101enum DataPart : unsigned {
102 DATA_PART_NONE = 0u, /**< Nothing */
103 DATA_PART_PROPERTIES = 1u, /**< Particle::p */
104 DATA_PART_POSITION = 2u, /**< Particle::r */
105 DATA_PART_MOMENTUM = 8u, /**< Particle::m */
106 DATA_PART_FORCE = 16u, /**< Particle::f */
107#ifdef ESPRESSO_BOND_CONSTRAINT
108 DATA_PART_RATTLE = 32u, /**< Particle::rattle */
109#endif
110 DATA_PART_BONDS = 64u /**< Particle::bonds */
112} // namespace Cells
113
114/**
115 * @brief Map the data parts flags from cells to those
116 * used internally by the ghost communication.
117 *
118 * @param data_parts data parts flags
119 * @return ghost communication flags
120 */
121unsigned map_data_parts(unsigned data_parts);
122
123namespace Cells {
124inline ParticleRange particles(std::span<Cell *const> cells) {
125 /* Find first non-empty cell */
126 auto first_non_empty = std::ranges::find_if(
127 cells, [](const Cell *c) { return not c->particles().empty(); });
128
129 return {CellParticleIterator(first_non_empty, cells.end()),
130 CellParticleIterator(cells.end())};
131}
132} // namespace Cells
133
134/**
135 * @brief Distance vector and length handed to pair kernels.
136 */
137struct Distance {
139 : vec21(vec21), dist2(vec21.norm2()) {}
140
142 double dist2;
143};
144
145namespace detail {
146// NOLINTNEXTLINE(bugprone-exception-escape)
147struct MinimalImageDistance {
148 BoxGeometry const box;
149
150 Distance operator()(Particle const &p1, Particle const &p2) const {
151 return Distance(box.get_mi_vector(p1.pos(), p2.pos()));
152 }
153};
154
155struct EuclidianDistance {
156 Distance operator()(Particle const &p1, Particle const &p2) const {
157 return Distance(p1.pos() - p2.pos());
158 }
159};
160} // namespace detail
161
162/** Describes a cell structure / cell system. Contains information
163 * about the communication of cell contents (particles, ghosts, ...)
164 * between different nodes and the relation between particle
165 * positions and the cell system. All other properties of the cell
166 * system which are not common between different cell systems have to
167 * be stored in separate structures.
168 */
169class CellStructure : public System::Leaf<CellStructure> {
170public:
171 static constexpr auto vector_length = 1;
172 struct AoSoA_pack;
175 using memory_space = Kokkos::HostSpace;
176 using ListAlgorithm = Cabana::HalfNeighborTag;
177 using ListType =
178 CustomVerletList<Kokkos::HostSpace, ListAlgorithm, Cabana::VerletLayout2D,
179 Cabana::TeamVectorOpTag>;
180
181private:
182 /** The local id-to-particle index */
183 std::vector<Particle *> m_particle_index;
184 /** Implementation of the primary particle decomposition */
185 std::unique_ptr<ParticleDecomposition> m_decomposition;
186 /** Active type in m_decomposition */
188 /** One of @ref Cells::Resort, announces the level of resort needed.
189 */
190 unsigned m_resort_particles = Cells::RESORT_NONE;
191 bool m_verlet_skin_set = false;
192 bool m_rebuild_verlet_list = true;
193 bool m_rebuild_verlet_list_cabana = true;
194 std::vector<std::pair<Particle *, Particle *>> m_verlet_list;
195 double m_le_pos_offset_at_last_resort = 0.;
196 /** @brief Verlet list skin. */
197 double m_verlet_skin = 0.;
198 double m_verlet_reuse = 0.;
199 int m_cached_max_local_particle_id = 0;
200 std::size_t m_num_local_particles_cached = 0;
201 int m_max_id = 0;
202 std::unique_ptr<Kokkos::View<int *>> m_id_to_index;
203 std::unique_ptr<ForceType> m_local_force;
204#ifdef ESPRESSO_ROTATION
205 std::unique_ptr<ForceType> m_local_torque;
206#endif
207#ifdef ESPRESSO_NPT
208 std::unique_ptr<VirialType> m_local_virial;
209#endif
210 std::unique_ptr<LocalBondState> m_bond_state;
211 std::unique_ptr<ListType> m_verlet_list_cabana;
212 /** particle properties using individual Kokkos Views */
213 std::unique_ptr<AoSoA_pack> m_aosoa;
214 /** The local id-to-index for aosoa data */
215 std::vector<Particle *> m_unique_particles;
216 std::shared_ptr<KokkosHandle> m_kokkos_handle;
217
218public:
219 CellStructure(BoxGeometry const &box);
220 virtual ~CellStructure();
221
222 bool use_verlet_list = true;
223
224 /**
225 * @brief Update local particle index.
226 *
227 * Update the entry for a particle in the local particle
228 * index.
229 *
230 * @param id Entry to update.
231 * @param p Pointer to the particle.
232 */
234 assert(id >= 0);
235 // cppcheck-suppress assertWithSideEffect
236 assert(not p or p->id() == id);
237
238 if (static_cast<unsigned int>(id) >= m_particle_index.size())
239 m_particle_index.resize(static_cast<unsigned int>(id + 1));
240
241 m_particle_index[static_cast<unsigned int>(id)] = p;
242 }
243
244 /**
245 * @brief Update local particle index.
246 *
247 * Update the entry for a particle in the local particle
248 * index.
249 *
250 * @param p Pointer to the particle.
251 */
253 update_particle_index(p.id(), std::addressof(p));
254 }
255
256 /**
257 * @brief Update local particle index.
258 *
259 * @param pl List of particles whose index entries should be updated.
260 */
262 for (auto &p : pl) {
263 update_particle_index(p.id(), std::addressof(p));
264 }
265 }
266
267 /**
268 * @brief Clear the particles index.
269 */
270 void clear_particle_index() { m_particle_index.clear(); }
271
272private:
273 /**
274 * @brief Append a particle to a list and update this
275 * particle index accordingly.
276 * @param pl List to add the particle to.
277 * @param p Particle to add.
278 */
279 Particle &append_indexed_particle(ParticleList &pl, Particle &&p) {
280 /* Check if cell may reallocate, in which case the index
281 * entries for all particles in this cell have to be
282 * updated. */
283 auto const may_reallocate = pl.size() >= pl.capacity();
284 auto &new_part = pl.insert(std::move(p));
285
286 if (may_reallocate)
288 else {
289 update_particle_index(new_part);
290 }
291
292 return new_part;
293 }
294
295public:
296 /**
297 * @brief Get a local particle by id.
298 *
299 * @param id Particle to get.
300 * @return Pointer to particle if it is local,
301 * nullptr otherwise.
302 */
304 assert(id >= 0);
305
306 if (static_cast<unsigned int>(id) >= m_particle_index.size())
307 return nullptr;
308
309 return m_particle_index[static_cast<unsigned int>(id)];
310 }
311
312 /** @overload */
313 const Particle *get_local_particle(int id) const {
314 assert(id >= 0);
315
316 if (static_cast<unsigned int>(id) >= m_particle_index.size())
317 return nullptr;
318
319 return m_particle_index[static_cast<unsigned int>(id)];
320 }
321
322 template <class InputRange, class OutputIterator>
323 void get_local_particles(InputRange ids, OutputIterator out) {
324 std::ranges::transform(ids, out,
325 [this](int id) { return get_local_particle(id); });
326 }
327
328 CellStructureType decomposition_type() const { return m_type; }
329
330 /** Maximal cutoff supported by current cell system. */
332
333 /** Maximal pair range supported by current cell system. */
335
337 return Cells::particles(decomposition().local_cells());
338 }
339
341 return Cells::particles(decomposition().ghost_cells());
342 }
343
344 std::size_t count_local_particles() const {
345 std::size_t count = 0;
346 for (auto const &cell : m_decomposition->local_cells()) {
347 count += cell->particles().size();
348 }
349 return count;
350 }
351
352 /** @brief whether to use parallel version of @ref for_each_local_particle */
353 bool use_parallel_for_each_local_particle() const { return true; }
354
355 /**
356 * @brief Run a kernel on all local particles.
357 * The kernel is assumed to be thread-safe.
358 */
360 bool parallel = true) const {
361 if (parallel and use_parallel_for_each_local_particle()) {
362 parallel_for_each_particle_impl(decomposition().local_cells(), f);
363 return;
364 }
365 for (auto &p : local_particles()) {
366 f(p);
367 }
368 }
369
370 /**
371 * @brief Run a kernel on all ghost particles.
372 * The kernel is assumed to be thread-safe.
373 */
375 for (auto &p : ghost_particles()) {
376 f(p);
377 }
378 }
379
380private:
381 /** Cell system dependent function to find the right cell for a
382 * particle.
383 * \param p Particle.
384 * \return pointer to cell where to put the particle, nullptr
385 * if the particle does not belong on this node.
386 */
387 Cell *particle_to_cell(const Particle &p) {
389 }
390 Cell const *particle_to_cell(const Particle &p) const {
392 }
393
394 void parallel_for_each_particle_impl(std::span<Cell *const> cells,
395 ParticleUnaryOp &f) const;
396
397public:
398 /**
399 * @brief Add a particle.
400 *
401 * Moves a particle into the cell system. This adds
402 * a particle to the local node, irrespective of where
403 * it belongs.
404 *
405 * @param p Particle to add.
406 * @return Pointer to the particle in the cell
407 * system.
408 */
410
411 /**
412 * @brief Add a particle.
413 *
414 * Moves a particle into the cell system, if it
415 * belongs to this node. Otherwise this does not
416 * have an effect and the particle is discarded.
417 * This can be used to add a particle without
418 * knowledge where it should be placed by calling
419 * the function on all nodes, it will then add
420 * the particle in exactly one place.
421 *
422 * @param p Particle to add.
423 * @return Pointer to particle if it is local, null
424 * otherwise.
425 */
427
428 /**
429 * @brief Remove a particle.
430 *
431 * Removes a particle and all bonds pointing
432 * to it. This is a collective call.
433 *
434 * @param id Id of particle to remove.
435 */
436 void remove_particle(int id);
437
438 /**
439 * @brief Get the maximal particle id on this node.
440 *
441 * This returns the highest particle id on
442 * this node, or -1 if there are no particles on this node.
443 */
444 int get_max_local_particle_id() const;
446 return m_cached_max_local_particle_id;
447 }
448 std::size_t get_num_local_particles_cached() const {
449 return m_num_local_particles_cached;
450 }
451 int get_local_pair_bond_numbers() const;
454 void set_local_bond_numbers(int p, int a, int d);
455#ifdef ESPRESSO_COLLISION_DETECTION
456 void clear_new_bonds();
457 void add_new_bond(int bond_id, std::vector<int> const &particle_ids);
458 void rebuild_bond_list();
459#endif // ESPRESSO_COLLISION_DETECTION
460
461 /**
462 * @brief Remove all particles from the cell system.
463 *
464 * This allows linear time removal of all particles from
465 * the system, removing each particle individually would
466 * be quadratic.
467 */
469
470 /**
471 * @brief Get the underlying particle decomposition.
472 *
473 * Should be used solely for informative purposes.
474 *
475 * @return The active particle decomposition.
476 */
478 return assert(m_decomposition), *m_decomposition;
479 }
480
481private:
483 return assert(m_decomposition), *m_decomposition;
484 }
485
486public:
487 /**
488 * @brief Increase the local resort level at least to @p level.
489 */
491 m_resort_particles |= level;
492 assert(m_resort_particles >= level);
493 }
494
495 /**
496 * @brief Get the currently scheduled resort level.
497 */
498 unsigned get_resort_particles() const { return m_resort_particles; }
499
500 /**
501 * @brief Set the resort level to sorted.
502 */
503 void clear_resort_particles() { m_resort_particles = Cells::RESORT_NONE; }
504
505 /**
506 * @brief Check whether a particle has moved further than half the skin
507 * since the last Verlet list update, thus requiring a resort.
508 * @param additional_offset Offset which is added to the distance the
509 * particle has travelled when comparing to half
510 * the Verlet skin (e.g., for Lees-Edwards BC).
511 * @return Whether a resort is needed.
512 */
513 bool
514 check_resort_required(Utils::Vector3d const &additional_offset = {}) const;
515
517 return m_le_pos_offset_at_last_resort;
518 }
519
520 /**
521 * @brief Synchronize number of ghosts.
522 */
523 void ghosts_count();
524
525 /**
526 * @brief Update ghost particles.
527 *
528 * Update ghost particles with data from the real particles.
529 *
530 * @param data_parts Particle parts to update, combination of @ref
531 * Cells::DataPart
532 */
533 void ghosts_update(unsigned data_parts);
534
535 /**
536 * @brief Update ghost particles, with particle resort if needed.
537 *
538 * Update ghost particles with data from the real particles.
539 * Resort particles if a resort is due.
540 *
541 * @param data_parts Particle parts to update, combination of @ref
542 * Cells::DataPart
543 */
544 void update_ghosts_and_resort_particle(unsigned data_parts);
545
546 /**
547 * @brief Add forces and torques from ghost particles to real particles.
548 */
550
551 /** Set forces and torques on all ghosts to zero. */
554 }
555
556#ifdef ESPRESSO_BOND_CONSTRAINT
557 /**
558 * @brief Add rattle corrections from ghost particles to real particles.
559 */
561#endif
562
563 /**
564 * @brief Resort particles.
565 */
566 void resort_particles(bool global_flag);
567
568 /** @brief Whether the Verlet skin is set. */
569 auto is_verlet_skin_set() const { return m_verlet_skin_set; }
570
571 /** @brief Get the Verlet skin. */
572 auto get_verlet_skin() const { return m_verlet_skin; }
573
574 /** @brief Set the Verlet skin. */
575 void set_verlet_skin(double value);
576
577 /** @brief Set the Verlet skin using a heuristic. */
579
580 void update_verlet_stats(int n_steps, int n_verlet_updates) {
581 if (n_verlet_updates > 0) {
582 m_verlet_reuse = n_steps / static_cast<double>(n_verlet_updates);
583 } else {
584 m_verlet_reuse = 0.;
585 }
586 }
587
588 /** @brief Average number of integration steps the Verlet list was re-used */
589 auto get_verlet_reuse() const { return m_verlet_reuse; }
590
591 /**
592 * @brief Resolve ids to particles.
593 *
594 * @throws BondResolutionError if one of the ids
595 * was not found.
596 *
597 * @param partner_ids Ids to resolve.
598 * @return Vector of Particle pointers.
599 */
600 auto resolve_bond_partners(std::span<const int> partner_ids) {
601 boost::container::static_vector<Particle *, 4> partners;
602 get_local_particles(partner_ids, std::back_inserter(partners));
603
604 /* Check if id resolution failed for any partner */
605 if (std::ranges::find(partners, nullptr) != partners.end()) {
606 throw BondResolutionError{};
607 }
608
609 return partners;
610 }
611
612private:
613 /**
614 * @brief Execute kernel for every bond on particle.
615 * @tparam Handler Callable, which can be invoked with
616 * (Particle, int, std::span<Particle *>),
617 * returning a bool.
618 * @param p Particles for whom the bonds are evaluated.
619 * @param handler is called for every bond, and handed
620 * p, the bond id and a span with the bond
621 * partners as arguments. Its return value
622 * should indicate if the bond was broken.
623 */
624 template <class Handler>
625 void execute_bond_handler(Particle &p, Handler handler) {
626 for (const BondView bond : p.bonds()) {
627 auto const partner_ids = bond.partner_ids();
628
629 try {
630 auto partners = resolve_bond_partners(partner_ids);
631 auto const partners_span = std::span(partners.data(), partners.size());
632 auto const bond_broken = handler(p, bond.bond_id(), partners_span);
633 if (bond_broken) {
634 bond_broken_error(p.id(), partner_ids);
635 }
636 } catch (const BondResolutionError &) {
637 bond_resolution_error(partner_ids);
638 }
639 }
640 }
641
642 /**
643 * @brief Go through ghost cells and remove the ghost entries from the
644 * local particle index.
645 */
646 void invalidate_ghosts() {
647 for (auto const &p : ghost_particles()) {
648 if (get_local_particle(p.id()) == &p) {
649 update_particle_index(p.id(), nullptr);
650 }
651 }
652 }
653
654 /** @brief Set the particle decomposition, keeping the particles. */
655 void set_particle_decomposition(
656 std::unique_ptr<ParticleDecomposition> &&decomposition) {
658
659 /* Swap in new cell system */
660 std::swap(m_decomposition, decomposition);
661
662 /* Add particles to new system */
663 for (auto &p : Cells::particles(decomposition->local_cells())) {
664 add_particle(std::move(p));
665 }
666 }
667
668public:
669 /**
670 * @brief Set the particle decomposition to @ref AtomDecomposition.
671 */
673
674 /**
675 * @brief Set the particle decomposition to @ref RegularDecomposition.
676 *
677 * @param range Interaction range.
678 * @param fully_connected_boundary neighbor cell directions for Lees-Edwards.
679 */
681 double range,
682 std::optional<std::pair<int, int>> fully_connected_boundary);
683
684 /**
685 * @brief Set the particle decomposition to @ref HybridDecomposition.
686 *
687 * @param cutoff_regular Interaction cutoff_regular.
688 * @param n_square_types Particle types to put into n_square decomposition.
689 */
690 void set_hybrid_decomposition(double cutoff_regular,
691 std::set<int> n_square_types);
692
693private:
694 /**
695 * @brief Run link_cell algorithm for local cells.
696 *
697 * @tparam Kernel Needs to be callable with (Particle, Particle, Distance).
698 * @param kernel Pair kernel functor.
699 */
700 template <class Kernel> void link_cell(Kernel kernel) {
701 auto const maybe_box = decomposition().minimum_image_distance();
702 auto const local_cells_span = decomposition().local_cells();
703 auto const first = boost::make_indirect_iterator(local_cells_span.begin());
704 auto const last = boost::make_indirect_iterator(local_cells_span.end());
705
706 if (maybe_box) {
708 first, last,
709 [&kernel, df = detail::MinimalImageDistance{decomposition().box()}](
710 Particle &p1, Particle &p2) { kernel(p1, p2, df(p1, p2)); });
711 } else {
712 if (decomposition().box().type() != BoxType::CUBOID) {
713 throw std::runtime_error("Non-cuboid box type is not compatible with a "
714 "particle decomposition that relies on "
715 "EuclideanDistance for distance calculation.");
716 }
718 first, last,
719 [&kernel, df = detail::EuclidianDistance{}](
720 Particle &p1, Particle &p2) { kernel(p1, p2, df(p1, p2)); });
721 }
722 }
723
724public:
725 auto get_max_id() const { return m_max_id; }
726
727 void set_kokkos_handle(std::shared_ptr<KokkosHandle> handle);
728 void rebuild_local_properties(double pair_cutoff);
730 void reset_local_force();
731
732 auto &get_id_to_index() { return *m_id_to_index; }
733 auto &get_local_force() { return *m_local_force; }
734#ifdef ESPRESSO_ROTATION
735 auto &get_local_torque() { return *m_local_torque; }
736#endif
737#ifdef ESPRESSO_NPT
738 auto &get_local_virial() { return *m_local_virial; }
739#endif
740 auto &get_aosoa() { return *m_aosoa; }
741 auto const &get_unique_particles() const { return m_unique_particles; }
742 auto const &get_verlet_list_cabana() const { return *m_verlet_list_cabana; }
743 auto &bond_state() { return *m_bond_state; }
744 auto const &bond_state() const { return *m_bond_state; }
747
748 [[nodiscard]] auto is_verlet_list_cabana_rebuild_needed() const {
749 return m_rebuild_verlet_list_cabana;
750 }
751
752 /**
753 * @brief Update bond storage(m_*_bond_list_kokkos and m_*_bond_id_kokkos).
754 * @param pair_count Index for pair bond storage.
755 * @param angle_count Index for angle bond storage.
756 * @param dihedral_count Index for dihedral bond storage.
757 * @param p Particle pointer.
758 */
759 void update_bond_storage(int &pair_count, int &angle_count,
760 int &dihedral_count, Particle const &p);
761
762 /**
763 * @brief Reset local properties of the Verlet list.
764 * @param cutoff Pair interaction cutoff.
765 * @return True if a rebuild is needed.
766 */
767 [[nodiscard]] auto prepare_verlet_list_cabana(double cutoff) {
768 auto const rebuild = is_verlet_list_cabana_rebuild_needed();
769 if (rebuild) {
770 // If we have to rebuild, we need to count the particles
771 set_index_map(); // parallelized index_map
772 // Create essential variables for MD
774 } else {
775 // If we do not rebuild we can use the saved map
777 }
778 return rebuild;
779 }
780
781 void rebuild_verlet_list_cabana(auto &&kernel, bool rebuild_verlet_list) {
783 if (rebuild_verlet_list) {
784 kernel(m_decomposition->local_cells(), m_decomposition->box(),
785 *m_verlet_list_cabana);
786 }
787 m_rebuild_verlet_list_cabana = false;
788 }
789
790 void set_index_map();
791
792 inline void cell_list_loop(auto &&kernel) {
793 kernel(m_decomposition->local_cells(), m_decomposition->box());
794 }
795
796private:
797 /** Non-bonded pair loop with verlet lists.
798 *
799 * @param pair_kernel Kernel to apply
800 * @param verlet_criterion Filter for verlet lists.
801 */
802 template <class PairKernel, class VerletCriterion>
803 void verlet_list_loop(PairKernel pair_kernel,
804 const VerletCriterion &verlet_criterion) {
805 /* In this case the verlet list update is attached to
806 * the pair kernel, and the verlet list is rebuilt as
807 * we go. */
808 if (m_rebuild_verlet_list) {
809 m_verlet_list.clear();
810
811 link_cell([&](Particle &p1, Particle &p2, Distance const &d) {
812 if (verlet_criterion(p1, p2, d)) {
813 m_verlet_list.emplace_back(&p1, &p2);
814 pair_kernel(p1, p2, d);
815 }
816 });
817
818 m_rebuild_verlet_list = false;
819 m_rebuild_verlet_list_cabana = true;
820 } else {
821 auto const maybe_box = decomposition().minimum_image_distance();
822 /* In this case the pair kernel is just run over the verlet list. */
823 if (maybe_box) {
824 auto const distance_function =
825 detail::MinimalImageDistance{decomposition().box()};
826 for (auto const &[p1, p2] : m_verlet_list) {
827 pair_kernel(*p1, *p2, distance_function(*p1, *p2));
828 }
829 } else {
830 auto const distance_function = detail::EuclidianDistance{};
831 for (auto const &[p1, p2] : m_verlet_list) {
832 pair_kernel(*p1, *p2, distance_function(*p1, *p2));
833 }
834 }
835 }
836 }
837
838public:
839 /** Bonded pair loop.
840 * @param bond_kernel Kernel to apply
841 */
842 template <class BondKernel> void bond_loop(BondKernel const &bond_kernel) {
843 for (auto &p : local_particles()) {
844 execute_bond_handler(p, bond_kernel);
845 }
846 }
847
848 /** Non-bonded pair loop.
849 * @param pair_kernel Kernel to apply
850 */
851 template <class PairKernel> void non_bonded_loop(PairKernel pair_kernel) {
852 link_cell(pair_kernel);
853 }
854
855 /** Non-bonded pair loop with potential use
856 * of verlet lists.
857 * @param pair_kernel Kernel to apply
858 * @param verlet_criterion Filter for verlet lists.
859 */
860 template <class PairKernel, class VerletCriterion>
861 void non_bonded_loop(PairKernel pair_kernel,
862 const VerletCriterion &verlet_criterion) {
863 if (use_verlet_list) {
864 verlet_list_loop(pair_kernel, verlet_criterion);
865 } else {
866 /* No verlet lists, just run the kernel with pairs from the cells. */
867 link_cell(pair_kernel);
868 }
869 }
870
871 /**
872 * @brief Check that particle index is commensurate with particles.
873 *
874 * For each local particles is checked that has a correct entry
875 * in the particles index, and that there are no excess (non-existing)
876 * particles in the index.
877 */
878 void check_particle_index() const;
879
880 /**
881 * @brief Check that particles are in the correct cell.
882 *
883 * This checks for all local particles that the result
884 * of particles_to_cell is the cell the particles is
885 * actually in, e.g. that the particles are sorted according
886 * to particles_to_cell.
887 */
888 void check_particle_sorting() const;
889
890public:
891 /**
892 * @brief Find cell a particle is stored in.
893 *
894 * For local particles, this returns the cell they
895 * are stored in, otherwise nullptr is returned.
896 *
897 * @param p Particle to find cell for
898 * @return Cell for particle or nullptr.
899 */
901 assert(not get_resort_particles());
902
903 if (p.is_ghost()) {
904 return nullptr;
905 }
906
907 return particle_to_cell(p);
908 }
909
910 /**
911 * @brief Run kernel on all particles inside local cell and its neighbors.
912 *
913 * @param p Particle to find cell for
914 * @param kernel Function with signature <tt>double(Particle const&,
915 * Particle const&, Utils::Vector3d const&)</tt>
916 * @return false if cell is not found, otherwise true
917 */
918 template <class Kernel>
920 Kernel &kernel) {
921 auto const cell = find_current_cell(p);
922
923 if (cell == nullptr) {
924 return false;
925 }
926
927 auto const maybe_box = decomposition().minimum_image_distance();
928
929 if (maybe_box) {
930 auto const distance_function =
931 detail::MinimalImageDistance{decomposition().box()};
932 short_range_neighbor_loop(p, cell, kernel, distance_function);
933 } else {
934 auto const distance_function = detail::EuclidianDistance{};
935 short_range_neighbor_loop(p, cell, kernel, distance_function);
936 }
937 return true;
938 }
939
940private:
941 template <class Kernel, class DistanceFunc>
942 void short_range_neighbor_loop(Particle const &p1, Cell *const cell,
943 Kernel &kernel, DistanceFunc const &df) {
944 /* Iterate over particles inside cell */
945 for (auto const &p2 : cell->particles()) {
946 if (p1.id() != p2.id()) {
947 auto const vec = df(p1, p2).vec21;
948 kernel(p1, p2, vec);
949 }
950 }
951 /* Iterate over all neighbors */
952 for (auto const neighbor : cell->neighbors().all()) {
953 /* Iterate over particles in neighbors */
954 if (neighbor != cell) {
955 for (auto const &p2 : neighbor->particles()) {
956 auto const vec = df(p1, p2).vec21;
957 kernel(p1, p2, vec);
958 }
959 }
960 }
961 }
962};
ParticleIterator< std::span< Cell *const >::iterator > CellParticleIterator
CellStructureType
Cell structure topology.
@ NSQUARE
Atom decomposition (N-square).
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.
void bond_broken_error(int id, std::span< const int > partner_ids)
void bond_resolution_error(std::span< const int > partner_ids)
Immutable view on a bond.
Definition BondList.hpp:44
ESPRESSO_ATTR_ALWAYS_INLINE Utils::Vector3< T > get_mi_vector(Utils::Vector3< T > const &a, Utils::Vector3< T > const &b) const
Get the minimum-image vector between two coordinates.
BoxType type() const
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 update_particle_index(ParticleList &pl)
Update local particle index.
void check_particle_sorting() const
Check that particles are in the correct cell.
void rebuild_verlet_list_cabana(auto &&kernel, bool rebuild_verlet_list)
auto & get_id_to_index()
std::size_t count_local_particles() const
auto const & bond_state() const
virtual ~CellStructure()
int get_local_angle_bond_numbers() const
void clear_resort_particles()
Set the resort level to sorted.
Cell * find_current_cell(const Particle &p)
Find cell a particle is stored in.
auto get_max_id() const
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.
int get_local_pair_bond_numbers() const
void clear_particle_index()
Clear the particles index.
auto prepare_verlet_list_cabana(double cutoff)
Reset local properties of the Verlet list.
static constexpr auto vector_length
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.
auto get_le_pos_offset_at_last_resort() const
int get_local_dihedral_bond_numbers() const
Kokkos::HostSpace memory_space
int get_cached_max_local_particle_id() const
void get_local_particles(InputRange ids, OutputIterator out)
void update_verlet_stats(int n_steps, int n_verlet_updates)
auto & get_local_torque()
void ghosts_reset_forces()
Set forces and torques on all ghosts to zero.
auto & get_local_virial()
void update_particle_index(int id, Particle *p)
Update local particle index.
void set_local_bond_numbers(int p, int a, int d)
void ghosts_reduce_forces()
Add forces and torques from ghost particles to real particles.
auto const & get_unique_particles() const
unsigned get_resort_particles() const
Get the currently scheduled resort level.
auto get_verlet_reuse() const
Average number of integration steps the Verlet list was re-used.
void rebuild_local_properties(double pair_cutoff)
void for_each_local_particle(ParticleUnaryOp &&f, bool parallel=true) const
Run a kernel on all local particles.
void non_bonded_loop(PairKernel pair_kernel)
Non-bonded pair loop.
void for_each_ghost_particle(ParticleUnaryOp &&f) const
Run a kernel on all ghost particles.
Utils::Vector3d max_range() const
Maximal pair range supported by current cell system.
void add_new_bond(int bond_id, std::vector< int > const &particle_ids)
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,...
const Particle * get_local_particle(int id) const
This is an overloaded member function, provided for convenience. It differs from the above function o...
auto resolve_bond_partners(std::span< const int > partner_ids)
Resolve ids to particles.
void bond_loop(BondKernel const &bond_kernel)
Bonded pair loop.
void ghosts_count()
Synchronize number of ghosts.
void set_resort_particles(Cells::Resort level)
Increase the local resort level at least to level.
void cell_list_loop(auto &&kernel)
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.
Cabana::HalfNeighborTag ListAlgorithm
auto get_verlet_skin() const
Get the Verlet skin.
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.
auto const & get_verlet_list_cabana() const
bool run_on_particle_short_range_neighbors(Particle const &p, Kernel &kernel)
Run kernel on all particles inside local cell and its neighbors.
bool use_parallel_for_each_local_particle() const
whether to use parallel version of for_each_local_particle
void remove_all_particles()
Remove all particles from the cell system.
ParticleRange local_particles() const
void update_particle_index(Particle &p)
Update local particle index.
void ghosts_reduce_rattle_correction()
Add rattle corrections from ghost particles to real particles.
CellStructureType decomposition_type() const
auto is_verlet_list_cabana_rebuild_needed() const
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 update_bond_storage(int &pair_count, int &angle_count, int &dihedral_count, Particle const &p)
Update bond storage(m_*_bond_list_kokkos and m_*_bond_id_kokkos).
void clear_bond_properties()
void non_bonded_loop(PairKernel pair_kernel, const VerletCriterion &verlet_criterion)
Non-bonded pair loop with potential use of verlet lists.
void reset_local_properties()
Definition Cell.hpp:96
auto & particles()
Particles.
Definition Cell.hpp:103
A distributed particle decomposition.
virtual Utils::Vector3d max_cutoff() const =0
Maximum supported cutoff.
virtual std::span< Cell *const > local_cells() const =0
Get pointer to local cells.
virtual Cell * particle_to_cell(Particle const &p)=0
Determine which cell a particle id belongs to.
virtual Utils::Vector3d max_range() const =0
Range in which calculations are performed.
virtual std::optional< BoxGeometry > minimum_image_distance() const =0
Return the box geometry needed for distance calculation if minimum image convention should be used ne...
virtual BoxGeometry const & box() const =0
A range of particles.
Abstract class that represents a component of the system.
std::size_t capacity() const
Capacity of the container.
Definition Bag.hpp:104
T & insert(T const &v)
Insert an element into the container.
Definition Bag.hpp:147
std::size_t size() const
Number of elements in the container.
Definition Bag.hpp:90
Returns true if the particles are to be considered for short range interactions.
Ghost particles and particle exchange.
void link_cell(CellIterator first, CellIterator last, PairKernel &&pair_kernel)
Iterates over all particles in the cell range, and over all pairs within the cells and with their nei...
Definition link_cell.hpp:32
DataPart
Flags to select particle parts for communication.
@ DATA_PART_MOMENTUM
Particle::m.
@ DATA_PART_FORCE
Particle::f.
@ DATA_PART_PROPERTIES
Particle::p.
@ DATA_PART_BONDS
Particle::bonds.
@ DATA_PART_NONE
Nothing.
@ DATA_PART_RATTLE
Particle::rattle.
@ DATA_PART_POSITION
Particle::r.
ParticleRange particles(std::span< Cell *const > cells)
auto constexpr new_part
Exception indicating that a particle id could not be resolved.
Distance vector and length handed to pair kernels.
Utils::Vector3d vec21
Distance(Utils::Vector3d const &vec21)
Struct holding all information for one particle.
Definition Particle.hpp:435
auto const & force_and_torque() const
Definition Particle.hpp:477
bool is_ghost() const
Definition Particle.hpp:480
auto const & pos() const
Definition Particle.hpp:471
auto const & id() const
Definition Particle.hpp:454