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