Loading [MathJax]/extensions/TeX/AMSmath.js
ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages Concepts
AtomDecomposition.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 "cell_system/Cell.hpp"
27
28#include "BoxGeometry.hpp"
29#include "Particle.hpp"
30#include "ghosts.hpp"
31
32#include <utils/Vector.hpp>
33
34#include <boost/mpi/communicator.hpp>
35
36#include <optional>
37#include <span>
38#include <utility>
39#include <vector>
40
41/**
42 * @brief Atom decomposition cell system.
43 *
44 * This implements a distributed particle storage
45 * by just evenly distributing the particles over
46 * all part-taking nodes. Pairs are found by just
47 * considering all pairs independent of logical or
48 * physical location, it has therefore quadratic time
49 * complexity in the number of particles.
50 *
51 * For a more detailed discussion please see @cite plimpton95a.
52 */
54 boost::mpi::communicator m_comm;
55 std::vector<Cell> cells;
56
57 std::vector<Cell *> m_local_cells;
58 std::vector<Cell *> m_ghost_cells;
59
60 GhostCommunicator m_exchange_ghosts_comm;
61 GhostCommunicator m_collect_ghost_force_comm;
62
63 BoxGeometry const &m_box;
64
65public:
66 AtomDecomposition(BoxGeometry const &m_box);
67 AtomDecomposition(boost::mpi::communicator comm, BoxGeometry const &box_geo);
68
69 void resort(bool global_flag, std::vector<ParticleChange> &diff) override;
70
71 GhostCommunicator const &exchange_ghosts_comm() const override {
72 return m_exchange_ghosts_comm;
73 }
75 return m_collect_ghost_force_comm;
76 }
77
78 std::span<Cell *const> local_cells() const override { return m_local_cells; }
79 std::span<Cell *const> ghost_cells() const override { return m_ghost_cells; }
80
81 /* Getter needed for HybridDecomposition */
82 auto const &get_local_cells() const { return m_local_cells; }
83 auto const &get_ghost_cells() const { return m_ghost_cells; }
84
85 /**
86 * @brief Determine which cell a particle id belongs to.
87 *
88 * Since there is only one local cell this is trivial.
89 *
90 * @param p Particle to find cell for.
91 * @return Pointer to cell or nullptr if not local.
92 */
93 Cell *particle_to_cell(Particle const &p) override {
94 return id_to_cell(p.id());
95 }
96 Cell const *particle_to_cell(Particle const &p) const override {
97 return id_to_cell(p.id());
98 }
99
100 Utils::Vector3d max_cutoff() const override;
101 Utils::Vector3d max_range() const override;
102
103 /* Return true if minimum image convention is
104 * needed for distance calculation. */
105 std::optional<BoxGeometry> minimum_image_distance() const override {
106 return m_box;
107 }
108
109 BoxGeometry const &box() const override { return m_box; }
110
111private:
112 /**
113 * @brief Find cell for id.
114 * @param id to find cell for.
115 * @return Cell for id.
116 */
117 Cell *id_to_cell(int id) {
118 return has_id(id) ? std::addressof(local()) : nullptr;
119 }
120
121 Cell const *id_to_cell(int id) const {
122 return has_id(id) ? std::addressof(local()) : nullptr;
123 }
124
125 /**
126 * @brief Get the local cell.
127 */
128 Cell &local() { return cells.at(static_cast<unsigned int>(m_comm.rank())); }
129 Cell const &local() const {
130 return cells.at(static_cast<unsigned int>(m_comm.rank()));
131 }
132
133 void configure_neighbors();
134 GhostCommunicator prepare_comm();
135
136 /**
137 * @brief Setup ghost communicators.
138 */
139 void configure_comms();
140
141 /**
142 * @brief Mark local and ghost cells.
143 */
144 void mark_cells();
145
146 /**
147 * @brief Determine which rank owns a particle id.
148 */
149 int id_to_rank(int id) const { return id % m_comm.size(); }
150
151 /**
152 * @brief Determine if this rank owns a particle id.
153 */
154 bool has_id(int id) const { return id_to_rank(id) == m_comm.rank(); }
155};
Vector implementation and trait types for boost qvm interoperability.
Atom decomposition cell system.
GhostCommunicator const & collect_ghost_force_comm() const override
Cell const * particle_to_cell(Particle const &p) const override
Utils::Vector3d max_range() const override
Utils::Vector3d max_cutoff() const override
GhostCommunicator const & exchange_ghosts_comm() const override
auto const & get_local_cells() const
std::optional< BoxGeometry > minimum_image_distance() const override
Cell * particle_to_cell(Particle const &p) override
Determine which cell a particle id belongs to.
auto const & get_ghost_cells() const
void resort(bool global_flag, std::vector< ParticleChange > &diff) override
std::span< Cell *const > local_cells() const override
BoxGeometry const & box() const override
std::span< Cell *const > ghost_cells() const override
Definition Cell.hpp:96
A distributed particle decomposition.
Ghost particles and particle exchange.
Properties for a ghost communication.
Definition ghosts.hpp:157
Struct holding all information for one particle.
Definition Particle.hpp:395
auto const & id() const
Definition Particle.hpp:414