ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
Loading...
Searching...
No Matches
RegularDecomposition.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 "cell_system/Cell.hpp"
27
28#include "BoxGeometry.hpp"
29#include "LocalBox.hpp"
30#include "Particle.hpp"
31#include "ParticleList.hpp"
32#include "ghosts.hpp"
33
34#include <utils/Vector.hpp>
35
36#include <boost/mpi/communicator.hpp>
37
38#include <cstddef>
39#include <optional>
40#include <span>
41#include <utility>
42#include <vector>
43
44/**
45 * @brief Regular decomposition cell system.
46 *
47 * The domain of a node is split into a 3D cell grid with dimension
48 * @ref RegularDecomposition::cell_grid "cell_grid". Together with one ghost
49 * cell layer on each side the overall dimension of the ghost cell grid is
50 * @ref RegularDecomposition::ghost_cell_grid "ghost_cell_grid". The regular
51 * decomposition enables the use of the linked cell algorithm
52 * which is in turn used for setting up the Verlet list for the
53 * system. You can see a 2D graphical representation of the linked
54 * cell grid below.
55 *
56 * \image html linked_cells.gif "Linked cells structure"
57 *
58 * 2D representation of a linked cell grid:
59 * <tt>cell_grid = {4,4}, ghost_cell_grid = {6,6}</tt>
60 *
61 * Each cell has @f$ 3^D @f$ neighbor cells. Since we deal with pair forces,
62 * it is sufficient to calculate only half of the interactions (Newton's law:
63 * action = reaction). We have chosen the upper half e.g. all neighbor
64 * cells with a higher linear index (for cell 14 they are marked in light
65 * blue). Caution: This implementation needs double sided ghost
66 * communication! For single sided ghost communication one would need
67 * some ghost-ghost cell interaction as well, which we do not need!
68 */
70 /** Grid dimensions per node. */
72 /** Cell size. */
74 /** Offset in global grid */
76 /** linked cell grid with ghost frame. */
78 /** inverse @ref RegularDecomposition::cell_size "cell_size". */
80
81 boost::mpi::communicator m_comm;
84 std::optional<std::pair<int, int>> m_fully_connected_boundary = {};
85 std::vector<Cell> cells;
86 std::vector<Cell *> m_local_cells;
87 std::vector<Cell *> m_ghost_cells;
90
91public:
92 RegularDecomposition(boost::mpi::communicator comm, double range,
93 BoxGeometry const &box_geo, LocalBox const &local_geo,
94 std::optional<std::pair<int, int>> fully_connected);
95
96 GhostCommunicator const &exchange_ghosts_comm() const override {
98 }
101 }
102
103 std::span<Cell *const> local_cells() const override { return m_local_cells; }
104 std::span<Cell *const> ghost_cells() const override { return m_ghost_cells; }
105
106 /* Getter needed for HybridDecomposition */
107 auto const &get_local_cells() const { return m_local_cells; }
108 auto const &get_ghost_cells() const { return m_ghost_cells; }
109
110 Cell *particle_to_cell(Particle const &p) override {
111 return position_to_cell(p.pos());
112 }
113
114 Cell const *particle_to_cell(Particle const &p) const override {
115 return position_to_cell(p.pos());
116 }
117
118 void resort(bool global, std::vector<ParticleChange> &diff) override;
119 Utils::Vector3d max_cutoff() const override;
120 Utils::Vector3d max_range() const override;
121
123
124 std::optional<BoxGeometry> minimum_image_distance() const override {
125 return {m_box};
126 }
127
128 BoxGeometry const &box() const override { return m_box; }
129
130private:
131 /** Fill @c m_local_cells list and @c m_ghost_cells list for use with regular
132 * decomposition.
133 */
134 void mark_cells();
135
136 /** Fill a communication cell pointer list. Fill the cell pointers of
137 * all cells which are inside a rectangular subgrid of the 3D cell
138 * grid starting from the
139 * lower left corner @p lc up to the high top corner @p hc. The cell
140 * pointer list @p part_lists must already be large enough.
141 * \param part_lists List of cell pointers to store the result.
142 * \param lc lower left corner of the subgrid.
143 * \param hc high up corner of the subgrid.
144 */
145 void fill_comm_cell_lists(ParticleList **part_lists,
146 Utils::Vector3i const &lc,
147 Utils::Vector3i const &hc);
148
149 int calc_processor_min_num_cells() const;
150
151 int position_to_cell_index(Utils::Vector3d const &pos) const;
152
153 /**
154 * @brief Get pointer to the cell which corresponds to the position if the
155 * position is in the node's spatial domain, otherwise a nullptr.
156 */
157 Cell *position_to_cell(Utils::Vector3d const &pos) {
158 auto const index = position_to_cell_index(pos);
159 return (index < 0) ? nullptr : &(cells.at(static_cast<std::size_t>(index)));
160 }
161 Cell const *position_to_cell(Utils::Vector3d const &pos) const {
162 auto const index = position_to_cell_index(pos);
163 return (index < 0) ? nullptr : &(cells.at(static_cast<std::size_t>(index)));
164 }
165
166 /**
167 * @brief Move particles into the cell system if it belongs to this node.
168 *
169 * Moves all particles from src into the local cell
170 * system if they do belong here. Otherwise the
171 * particles are moved into rest.
172 *
173 * @param src Particles to move.
174 * @param rest Output list for left-over particles.
175 * @param modified_cells Local cells that were touched.
176 */
177 void move_if_local(ParticleList &src, ParticleList &rest,
178 std::vector<ParticleChange> &modified_cells);
179
180 /**
181 * @brief Split particle list by direction.
182 *
183 * Moves all particles from @p src into @p left
184 * or @p right depending on whether they belong
185 * to the left or right side of the local node
186 * in direction @p dir.
187 *
188 * @param src Particles to sort.
189 * @param left Particles that should go to the left
190 * @param right Particles that should go to the right
191 * @param dir Direction to consider.
192 */
193 void move_left_or_right(ParticleList &src, ParticleList &left,
194 ParticleList &right, int dir) const;
195
196 /**
197 * @brief One round of particle exchange with the next neighbors.
198 *
199 * @param[in] pl Particle on the move
200 * @param[out] modified_cells Cells that got touched.
201 */
202 void exchange_neighbors(ParticleList &pl,
203 std::vector<ParticleChange> &modified_cells);
204
205 /**
206 * @brief Calculate cell grid dimensions, cell sizes and number of cells.
207 *
208 * Calculates the cell grid, based on the local box size and the range.
209 * If the number of cells is larger than @c max_num_cells,
210 * it increases @c max_range until the number of cells is
211 * smaller or equal to @c max_num_cells. It sets:
212 * @c cell_grid,
213 * @c ghost_cell_grid,
214 * @c cell_size, and
215 * @c inv_cell_size.
216 *
217 * @param range interaction range. All pairs closer
218 * than this distance are found.
219 */
220 void create_cell_grid(double range);
221
222 /** Init cell interactions for cell system regular decomposition.
223 * Initializes the interacting neighbor cell list of a cell.
224 * This list of interacting neighbor cells is used by the Verlet
225 * algorithm.
226 */
227 void init_cell_interactions();
228
229 /** Create communicators for cell structure regular decomposition (see \ref
230 * GhostCommunicator).
231 */
232 GhostCommunicator prepare_comm();
233
234 /** Maximal number of cells per node. In order to avoid memory
235 * problems due to the cell grid, one has to specify the maximal
236 * number of cells. If the number of cells is larger
237 * than @c max_num_cells, the cell grid is reduced.
238 * @c max_num_cells has to be larger than 27, e.g. one inner cell.
239 */
240 static constexpr int max_num_cells = 32768;
241};
Vector implementation and trait types for boost qvm interoperability.
Definition Cell.hpp:96
A distributed particle decomposition.
Ghost particles and particle exchange.
Properties for a ghost communication.
Definition ghosts.hpp:159
Struct holding all information for one particle.
Definition Particle.hpp:435
auto const & pos() const
Definition Particle.hpp:471
Regular decomposition cell system.
GhostCommunicator m_exchange_ghosts_comm
Utils::Vector3i ghost_cell_grid
linked cell grid with ghost frame.
std::optional< std::pair< int, int > > m_fully_connected_boundary
Utils::Vector3d max_cutoff() const override
Utils::Vector3d inv_cell_size
inverse cell_size.
std::vector< Cell * > m_ghost_cells
BoxGeometry const & m_box
Cell * particle_to_cell(Particle const &p) override
void resort(bool global, std::vector< ParticleChange > &diff) override
std::vector< Cell > cells
Utils::Vector3d cell_size
Cell size.
GhostCommunicator const & collect_ghost_force_comm() const override
std::span< Cell *const > local_cells() const override
std::vector< Cell * > m_local_cells
Cell const * particle_to_cell(Particle const &p) const override
BoxGeometry const & box() const override
auto const & get_local_cells() const
Utils::Vector3d max_range() const override
std::span< Cell *const > ghost_cells() const override
GhostCommunicator const & exchange_ghosts_comm() const override
GhostCommunicator m_collect_ghost_force_comm
Utils::Vector3i cell_grid
Grid dimensions per node.
auto const & get_ghost_cells() const
Utils::Vector3i cell_offset
Offset in global grid.
std::optional< BoxGeometry > minimum_image_distance() const override
boost::mpi::communicator m_comm