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