Utopia 2
Framework for studying models of complex & adaptive systems.
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Static Public Attributes | Protected Member Functions | Protected Attributes | Private Member Functions | Private Attributes | List of all members
Utopia::SquareGrid< Space > Class Template Reference

A grid discretization using square cells. More...

#include <square.hh>

Inheritance diagram for Utopia::SquareGrid< Space >:
Inheritance graph
[legend]
Collaboration diagram for Utopia::SquareGrid< Space >:
Collaboration graph
[legend]

Public Types

using Base = Grid< Space >
 Base class type.
 
using SpaceVec = typename Space::SpaceVec
 The type of vectors that have a relation to physical space.
 
using MultiIndex = MultiIndexType< dim >
 The type of multi-index like arrays, e.g. the grid shape.
 
using Config = DataIO::Config
 The configuration type.
 
- Public Types inherited from Utopia::Grid< Space >
using Self = Grid< Space >
 Type of this class, i.e. the base grid class.
 
using SpaceVec = typename Space::SpaceVec
 The type of vectors that have a relation to physical space.
 
using MultiIndex = MultiIndexType< dim >
 The type of multi-index like arrays, e.g. the grid shape.
 
using Config = DataIO::Config
 The configuration type.
 

Public Member Functions

 SquareGrid (std::shared_ptr< Space > space, const Config &cfg)
 Construct a rectangular grid discretization.
 
IndexType num_cells () const override
 Number of square cells required to fill the physical space.
 
SpaceVec effective_resolution () const override
 The effective cell resolution into each physical space dimension.
 
MultiIndex shape () const override
 Get shape of the square grid.
 
GridStructure structure () const override
 Structure of the grid.
 
MultiIndex midx_of (const IndexType id) const override
 Returns the multi-index of the cell with the given ID.
 
SpaceVec barycenter_of (const IndexType id) const override
 Returns the barycenter of the cell with the given ID.
 
SpaceVec extent_of (const IndexType) const override
 Returns the extent of the cell with the given ID.
 
std::vector< SpaceVecvertices_of (const IndexType id) const override
 Returns the vertices of the cell with the given ID.
 
IndexType cell_at (const SpaceVec &pos) const override
 Return the ID of the cell covering the given point in physical space.
 
std::set< IndexTypeboundary_cells (std::string select="all") const override
 Retrieve a set of cell indices that are at a specified boundary.
 
- Public Member Functions inherited from Utopia::Grid< Space >
 Grid (std::shared_ptr< Space > space, const Config &cfg)
 Construct a grid discretization.
 
virtual ~Grid ()=default
 Virtual destructor to allow polymorphic destruction.
 
IndexContainer neighbors_of (const IndexType id) const
 Returns the indices of the neighbors of the cell with the given ID.
 
void select_neighborhood (NBMode nb_mode, const Config &nb_params={})
 
const NBModenb_mode () const
 Const reference to the currently selected neighborhood mode.
 
const Confignb_params () const
 The neighborhood parameters of the currently selected neighborhood.
 
auto nb_size () const
 Maximum size of the currently selected neighborhood.
 
auto resolution () const
 Get scalar resolution value of this grid.
 
std::string structure_name () const
 Structure of the grid as std::string.
 
const std::shared_ptr< Space > & space () const
 Const reference to the space this grid maps to.
 
bool is_periodic () const
 Whether the space this grid maps to is periodic.
 

Static Public Attributes

static constexpr DimType dim = Space::dim
 The dimensionality of the space to be discretized (for easier access)
 
- Static Public Attributes inherited from Utopia::Grid< Space >
static constexpr DimType dim = Space::dim
 The dimensionality of the space to be discretized (for easier access)
 

Protected Member Functions

NBFuncID< Baseget_nb_func (NBMode nb_mode, const Config &nb_params) override
 Retrieve the neighborhood function depending on the mode and parameters.
 
NBFuncID< Baseget_nb_func_vonNeumann (const Config &nb_params)
 Returns a standalone von-Neumann neighborhood function.
 
NBFuncID< Baseget_nb_func_Moore (const Config &nb_params)
 Returns a standalone Moore neighborhood function.
 
template<DimType axis>
constexpr IndexType id_shift () const
 Return the shift in cell indices necessary if moving along an axis.
 
template<DimType axis, bool periodic>
void add_neighbors_in_ (const IndexType root_id, IndexContainer &neighbor_ids) const
 Add both direct neighbors to a container of indices.
 
template<DimType axis, bool periodic>
void add_low_val_neighbor_in_ (const IndexType root_id, const DistType distance, IndexContainer &neighbor_ids) const
 Add a neighbor on the low (ID) value side to an index container.
 
template<DimType axis, bool periodic>
void add_high_val_neighbor_in_ (const IndexType root_id, const DistType distance, IndexContainer &neighbor_ids) const
 Add a neighbor on the high (ID) value side to an index container.
 
DistType expected_num_neighbors (const NBMode &nb_mode, const Config &nb_params) const override
 Computes the expected number of neighbors for a neighborhood mode.
 
template<bool check_shape = false>
DistType get_nb_param_distance (const Config &params) const
 Extract the distance neighborhood parameter from the given config.
 

Protected Attributes

NBFuncID< Base_nb_vonNeumann_periodic
 The Von-Neumann neighborhood for periodic grids.
 
NBFuncID< Base_nb_vonNeumann_nonperiodic
 The Von-Neumann neighborhood for non-periodic grids.
 
NBFuncID< Base_nb_Moore_periodic
 Moore neighbors for periodic 2D grid.
 
NBFuncID< Base_nb_Moore_nonperiodic
 Moore neighbors for non-periodic 2D grid.
 
- Protected Attributes inherited from Utopia::Grid< Space >
const std::shared_ptr< Space_space
 The space that is to be discretized.
 
const DistType _resolution
 How many cells to place per length unit of space.
 
NBMode _nb_mode
 Neighborhood mode.
 
Config _nb_params
 Neighborhood parameters.
 
NBFuncID< Self_nb_func
 Neighborhood function (working on cell IDs)
 
NBFuncID< Self_nb_empty
 A neighborhood function for empty neighborhood.
 

Private Member Functions

MultiIndex determine_shape () const
 Given the resolution, return the grid shape required to fill the space.
 

Private Attributes

const MultiIndex _shape
 The (multi-index) shape of the grid, resulting from resolution.
 
const SpaceVec _cell_extent
 The extent of each cell of this square discretization (same for all)
 

Detailed Description

template<class Space>
class Utopia::SquareGrid< Space >

A grid discretization using square cells.

This is a grid discretization where the cells are vector spaces that are spanned by orthogonal basis vectors and each cell has the same physical extent in each dimension. In the 2D case, this refers to perfectly square cells; in 3D these would be perfect cubes, etc.

Note
Indices are constructed in "Fortran"-style, i.e. with the first index changing fastest (also called "column major"). For example, when iterating over the cell IDs and having a 2D space, the iteration goes first along the x-axis and then along the y-axis: 0, 1, …, N_x, N_x + 1, …, 2*N_x, 2*N_x + 1, …, N_x * N_y - 1

Member Typedef Documentation

◆ Base

Base class type.

◆ Config

The configuration type.

◆ MultiIndex

template<class Space >
using Utopia::SquareGrid< Space >::MultiIndex = MultiIndexType<dim>

The type of multi-index like arrays, e.g. the grid shape.

◆ SpaceVec

The type of vectors that have a relation to physical space.

Constructor & Destructor Documentation

◆ SquareGrid()

template<class Space >
Utopia::SquareGrid< Space >::SquareGrid ( std::shared_ptr< Space space,
const Config cfg 
)
inline

Construct a rectangular grid discretization.

Parameters
spaceThe space to construct the discretization for
cfgFurther configuration parameters
65 :
66 Base(space, cfg),
69 {
70 if constexpr (dim > 1) {
71 // Make sure the cells really are square
73
74 if (eff_res.min() != eff_res.max()) {
75 // Format effective resolution to show it in error message
76 std::stringstream efr_ss;
77 efr_ss << eff_res;
78
79 throw std::invalid_argument("Given the extent of the physical "
80 "space and the specified resolution, a mapping with "
81 "exactly square cells could not be found! Either adjust "
82 "the physical space, the resolution of the grid, or "
83 "choose another grid. Effective resolution was:\n"
84 + efr_ss.str() + ", but should be the same in all "
85 "dimensions!");
86 }
87 }
88 }
const std::shared_ptr< Space > & space() const
Const reference to the space this grid maps to.
Definition base.hh:299
const SpaceVec _cell_extent
The extent of each cell of this square discretization (same for all)
Definition square.hh:55
Grid< Space > Base
Base class type.
Definition square.hh:34
SpaceVec effective_resolution() const override
The effective cell resolution into each physical space dimension.
Definition square.hh:105
static constexpr DimType dim
The dimensionality of the space to be discretized (for easier access)
Definition square.hh:37
typename Space::SpaceVec SpaceVec
The type of vectors that have a relation to physical space.
Definition square.hh:40
MultiIndex determine_shape() const
Given the resolution, return the grid shape required to fill the space.
Definition square.hh:388
const MultiIndex _shape
The (multi-index) shape of the grid, resulting from resolution.
Definition square.hh:52
Container select_entities(const Manager &mngr, const DataIO::Config &sel_cfg)
Select entities according to parameters specified in a configuration.
Definition select.hh:213

Member Function Documentation

◆ add_high_val_neighbor_in_()

template<class Space >
template<DimType axis, bool periodic>
void Utopia::SquareGrid< Space >::add_high_val_neighbor_in_ ( const IndexType  root_id,
const DistType  distance,
IndexContainer neighbor_ids 
) const
inlineprotected

Add a neighbor on the high (ID) value side to an index container.

This function takes an index container and populates it with the index of a neighboring cell in different dimensions, specified by template parameter 0 < axis < number of dimensions - 1.

The algorithm first calculates whether the given root cell index has a back boundary in the chosen dimension. If so, the neighboring cell is only added if the grid is periodic.

Parameters
root_idWhich cell to find the agents of
distanceWhich distance the neighbor has to the root cell
neighbor_idsThe container to populate with the indices
Template Parameters
axisThe axis along which to add the neighbor (0-based!)
periodicWhether the grid is periodic
Returns
void
1035 {
1036 // Assure the number of dimensions is supported
1037 static_assert(dim <= 2,
1038 "Unsupported dimensionality of underlying space! Need be 1 or 2.");
1039 static_assert(axis < dim);
1040
1041 // If the distance is zero, no neighbor can be added; return nothing.
1042 if (distance == 0) {
1043 return;
1044 }
1045
1046 // Check if the neighbor to be added would pass a high value boundary.
1047 // Do so by computing a "normalized" ID along the desired dimension in
1048 // which the neighbor is to be added (always in [0, shape[d] - 1]) and
1049 // then compare to the distance from the high value boundary.
1051 >= _shape[axis] - distance)
1052 {
1053 if constexpr (periodic) {
1054 neighbor_ids.push_back( root_id
1055 + distance * id_shift<axis>()
1056 - id_shift<axis+1>());
1057 }
1058 }
1059 else {
1060 neighbor_ids.push_back( root_id
1061 + distance * id_shift<axis>());
1062 }
1063 }

◆ add_low_val_neighbor_in_()

template<class Space >
template<DimType axis, bool periodic>
void Utopia::SquareGrid< Space >::add_low_val_neighbor_in_ ( const IndexType  root_id,
const DistType  distance,
IndexContainer neighbor_ids 
) const
inlineprotected

Add a neighbor on the low (ID) value side to an index container.

This function takes an index container and populates it with the indices of neighboring cells in different dimensions, specified by template parameter 0 < axis < number of dimensions - 1.

The algorithm first calculates whether the given root cell index has a front boundary in the chosen dimension. If so, the neighboring cell is only added if the grid is periodic.

Parameters
root_idWhich cell to find the agents of
distanceWhich distance the neighbor has to the root cell
neighbor_idsThe container to populate with the indices
Template Parameters
axisThe axis along which to add the neighbor (0-based!)
periodicWhether the grid is periodic
Returns
void
983 {
984 // Assure the number of dimensions is supported
985 static_assert(dim <= 2,
986 "Unsupported dimensionality of underlying space! Need be 1 or 2.");
987 static_assert(axis < dim);
988
989 // If the distance is zero, no neighbor can be added; return nothing.
990 if (distance == 0) {
991 return;
992 }
993
994 // Check if the neighbor to be added would pass a low value boundary.
995 // Do so by computing a "normalized" ID along the desired dimension in
996 // which the neighbor is to be added (always in [0, shape[d] - 1]) and
997 // then compare to the distance
999 < distance)
1000 {
1001 if constexpr (periodic) {
1002 neighbor_ids.push_back( root_id
1003 - distance * id_shift<axis>()
1004 + id_shift<axis+1>());
1005 }
1006 }
1007 else {
1008 neighbor_ids.push_back( root_id
1009 - distance * id_shift<axis>());
1010 }
1011 }

◆ add_neighbors_in_()

template<class Space >
template<DimType axis, bool periodic>
void Utopia::SquareGrid< Space >::add_neighbors_in_ ( const IndexType  root_id,
IndexContainer neighbor_ids 
) const
inlineprotected

Add both direct neighbors to a container of indices.

This function takes an index container and populates it with the indices of neighboring cells in different dimensions, specified by template parameter 0 < axis < number of dimensions - 1.

The algorithm first calculates whether the given root cell index has a front or back boundary in the chosen dimension. If so, the neighboring cell is only added if the grid is periodic.

Parameters
root_idWhich cell to find the agents of
neighbor_idsThe container to populate with the indices
Template Parameters
axisThe axis along which to add the neighbors (0-based!)
periodicWhether the grid is periodic
Returns
void
922 {
923 // Assure the number of dimensions is supported
924 static_assert(dim <= 2,
925 "Unsupported dimensionality of underlying space! Need be 1 or 2.");
926 static_assert(axis < dim);
927
928 // Compute a "normalized" ID along the desired dimension in which the
929 // neighbors are to be added. Is always in [0, shape[d] - 1].
930 const auto nrm_id = ( (root_id % id_shift<axis+1>())
931 / id_shift<axis>());
932 // NOTE _Should_ also work for d > 2, but need tests for that.
933
934 // Check if at low value boundary
935 if (nrm_id == 0) {
936 if constexpr (periodic) {
937 neighbor_ids.push_back( root_id
939 + id_shift<axis+1>());
940 }
941 // else: not periodic; nothing to add here
942 }
943 else {
944 // Not at boundary; no need for the correction term
945 neighbor_ids.push_back(root_id - id_shift<axis>());
946 }
947
948 // Check if at high value boundary
949 if (nrm_id == _shape[axis] - 1) {
950 if constexpr (periodic) {
951 neighbor_ids.push_back( root_id
953 - id_shift<axis+1>());
954 }
955 }
956 else {
957 neighbor_ids.push_back(root_id + id_shift<axis>());
958 }
959 }

◆ barycenter_of()

template<class Space >
SpaceVec Utopia::SquareGrid< Space >::barycenter_of ( const IndexType  id) const
inlineoverridevirtual

Returns the barycenter of the cell with the given ID.

Note
This method does not perform bounds checking of the given ID!

Implements Utopia::Grid< Space >.

140 {
141 // Offset on grid + offset within cell
142 return (midx_of(id) % _cell_extent) + (_cell_extent/2.);
143 // NOTE The %-operator performs element-wise multiplication
144 }
MultiIndex midx_of(const IndexType id) const override
Returns the multi-index of the cell with the given ID.
Definition square.hh:125

◆ boundary_cells()

template<class Space >
std::set< IndexType > Utopia::SquareGrid< Space >::boundary_cells ( std::string  select = "all") const
inlineoverridevirtual

Retrieve a set of cell indices that are at a specified boundary.

Note
For a periodic space, an empty container is returned; no error or warning is emitted.
Parameters
selectWhich boundary to return the cell IDs of. If 'all', all boundary cells are returned. Other available values depend on the dimensionality of the grid: 1D: left, right 2D: bottom, top 3D: back, front

Implements Utopia::Grid< Space >.

258 {
259 static_assert(dim <= 2,
260 "SquareGrid::boundary_cells only implemented for 1D and 2D!");
261
262 // else: non-periodic space
263
264 // The target set all IDs are to be emplaced in
265 std::set<IndexType> bc_ids;
266
267 // Distinguish by dimensionality of the space
268 if constexpr (dim == 1) {
269 if (select != "all" and select != "left" and select != "right") {
270 throw std::invalid_argument("Invalid value for argument "
271 "`select` in call to method SquareGrid::boundary_cells! "
272 "Available arguments (for currently selected "
273 "dimensionality) are: "
274 "'all', 'left', 'right'. Given value: '" + select + "'");
275 }
276 // For periodic space, this is easy:
277 if (this->is_periodic()) {
278 return {};
279 }
280
281 // Left boundary (consists only of one cell)
282 if (select == "all" or select == "left") {
283 bc_ids.emplace(0);
284 }
285
286 // Right boundary (also consists only of one cell)
287 if (select == "all" or select == "right") {
288 bc_ids.emplace_hint(bc_ids.end(), _shape[0] - 1);
289 }
290 }
291 else if constexpr (dim == 2) {
292 if ( select != "all"
293 and select != "left" and select != "right"
294 and select != "bottom" and select != "top")
295 {
296 throw std::invalid_argument("Invalid value for argument "
297 "`select` in call to method SquareGrid::boundary_cells! "
298 "Available arguments (for currently selected "
299 "dimensionality) are: "
300 "'all', 'left', 'right', 'bottom', 'top'. Given value: '"
301 + select + "'");
302 }
303
304 // For periodic space, this is easy:
305 if (this->is_periodic()) {
306 return {};
307 }
308
309 // NOTE It is important to use the hinting features of std::set
310 // here, which allow to run the following in amortized
311 // constant time instead of logarithmic with set size.
312 // Below, it always makes sense to hint at inserting right
313 // before the end.
314 // Hint for the first element needs to be the beginning
315 auto hint = bc_ids.begin();
316
317 // Bottom boundary (lowest IDs)
318 if (select == "all" or select == "bottom") {
319 // 0, ..., _shape[0] - 1
320 for (DistType id = 0; id < _shape[0]; id++) {
321 bc_ids.emplace_hint(hint, id);
322 hint = bc_ids.end();
323 }
324 }
325
326 // Left boundary
327 if (select == "left") {
328 // First IDs in _shape[1] rows: 0, _shape[0], 2*_shape[0], ...
329 for (DistType row = 0; row < _shape[1]; row++) {
330 bc_ids.emplace_hint(hint, row * _shape[0]);
331 hint = bc_ids.end();
332 }
333 }
334
335 // Right boundary
336 if (select == "right") {
337 // Last IDs in _shape[1] rows
338 const auto offset = _shape[0] - 1;
339
340 for (DistType row = 0; row < _shape[1]; row++) {
341 bc_ids.emplace_hint(hint, offset + row * _shape[0]);
342 hint = bc_ids.end();
343 }
344 }
345
346 // Left AND right (only for 'all' case, allows better hints)
347 if (select == "all") {
348 // First and last IDs in _shape[1] rows
349 const auto offset = _shape[0] - 1;
350
351 for (DistType row = 0; row < _shape[1]; row++) {
352 // Left boundary cell
353 bc_ids.emplace_hint(hint, row * _shape[0]);
354
355 // Right boundary cell (higher than left cell ID)
356 bc_ids.emplace_hint(bc_ids.end(),
357 offset + row * _shape[0]);
358
359 hint = bc_ids.end();
360 }
361 }
362
363 // Top boundary (highest IDs)
364 if (select == "all" or select == "top") {
365 // _shape[0] * (_shape[1]-1), ..., _shape[0] * _shape[1] - 1
366 for (DistType id = _shape[0] * (_shape[1]-1);
367 id < _shape[0] * _shape[1]; id++)
368 {
369 bc_ids.emplace_hint(hint, id);
370 hint = bc_ids.end();
371 }
372 }
373 }
374
375 return bc_ids;
376 }
bool is_periodic() const
Whether the space this grid maps to is periodic.
Definition base.hh:304
unsigned int DistType
Type for distancens, i.e. intermediately long unsigned integers.
Definition types.hh:37

◆ cell_at()

template<class Space >
IndexType Utopia::SquareGrid< Space >::cell_at ( const SpaceVec pos) const
inlineoverridevirtual

Return the ID of the cell covering the given point in physical space.

Cells are interpreted as covering half-open intervals in space, i.e., including their low-value edges and excluding their high-value edges. The special case of points on high-value edges for non-periodic space behaves such that these points are associated with the cells at the boundary.

Note
This function always returns IDs of cells that are inside physical space. For non-periodic space, a check is performed whether the given point is inside the physical space associated with this grid. For periodic space, the given position is mapped back into the physical space.

Implements Utopia::Grid< Space >.

190 {
191 static_assert(dim == 2,
192 "SquareGrid::cell_at only implemented for 2D!");
193
194 // The multi-index element type to use for static casts
195 using midx_et = typename MultiIndex::elem_type;
196
197 // The multi-index to be calculated
199
200 // Distinguish periodic and non-periodic case
201 // NOTE While there is some duplicate code, this is the configuration
202 // with the least amount of unnecessary / duplicate value checks.
203 if (this->is_periodic()) {
204 // Calculate the real-valued position in units of cell extents,
205 // using the position mapped back into space. That function takes
206 // care to map the high-value boundary to the low-value boundary.
207 const SpaceVec ridx = ( this->_space->map_into_space(pos)
208 / _cell_extent);
209
210 // Can now calculate the integer multi index, rounding down
211 midx = {static_cast<midx_et>(ridx[0]),
212 static_cast<midx_et>(ridx[1])};
213 }
214 else {
215 // Make sure the given coordinate is actually within the space
216 if (not this->_space->contains(pos)) {
217 throw std::invalid_argument("The given position is outside "
218 "the non-periodic space associated with this grid!");
219 }
220
221 // Calculate the real-valued position in units of cell extents
222 const SpaceVec ridx = pos / _cell_extent;
223
224 // Calculate the integer multi index, rounding down
225 midx = {static_cast<midx_et>(ridx[0]),
226 static_cast<midx_et>(ridx[1])};
227
228 // Associate points on high-value boundaries with boundary cells
229 if (midx[0] == _shape[0]) {
230 midx[0]--;
231 }
232 if (midx[1] == _shape[1]) {
233 midx[1]--;
234 }
235 // Note that with _shape having only elements >= 1, the decrement
236 // does not cause any issues.
237 }
238
239 // From the multi index, calculate the corresponding ID
240 return midx[0] + (midx[1] * _shape[0]);
241 // Equivalent to:
242 // midx[0] * id_shift<0>()
243 // + midx[1] * id_shift<1>()
244 }
const std::shared_ptr< Space > _space
The space that is to be discretized.
Definition base.hh:120
MultiIndexType< dim > MultiIndex
The type of multi-index like arrays, e.g. the grid shape.
Definition square.hh:43

◆ determine_shape()

template<class Space >
MultiIndex Utopia::SquareGrid< Space >::determine_shape ( ) const
inlineprivate

Given the resolution, return the grid shape required to fill the space.

Integer rounding takes place here. A physical space of extents of 2.1 length units in each dimension and a resolution of two cells per unit length will result in 4 cells in each dimension, each cell's size scaled up slightly and the effective resolution thus slightly smaller than the specified resolution.

388 {
390
391 for (DimType i = 0; i < dim; i++) {
392 shape[i] = this->_space->extent[i] * this->_resolution;
393 }
394
395 return shape;
396 }
const DistType _resolution
How many cells to place per length unit of space.
Definition base.hh:127
MultiIndex shape() const override
Get shape of the square grid.
Definition square.hh:111
unsigned short DimType
Type for dimensions, i.e. very small unsigned integers.
Definition types.hh:34

◆ effective_resolution()

template<class Space >
SpaceVec Utopia::SquareGrid< Space >::effective_resolution ( ) const
inlineoverridevirtual

The effective cell resolution into each physical space dimension.

For a square lattice, this is just the quotient of grid shape and extent of physical space, separately in each dimension

Implements Utopia::Grid< Space >.

105 {
106 // Use element-wise division by the physical extent (double)
107 return _shape / this->_space->extent;
108 }

◆ expected_num_neighbors()

template<class Space >
DistType Utopia::SquareGrid< Space >::expected_num_neighbors ( const NBMode nb_mode,
const Config nb_params 
) const
inlineoverrideprotectedvirtual

Computes the expected number of neighbors for a neighborhood mode.

This function is used to calculate the amount of memory that should be reserved for the neighbor_ids vector. For the calculation it uses the member variables: dim and the distance parameter from the given configuration.

For a von Neumann neighborhood, the number of neighbors is: { 2 * distance for distance = 1 N(dim, distance) = { 2 * sum_{distances} N(dim-1, distance) { for distance > 1)

For a Moore neighborhood, the number of neighbors is:

N(dim, distance) = (2 * distance + 1)^2 - 1

Warning
The expected number of neighbors is not equal to the actually calculated number of neighbors. For example in a nonperiodic grids the expected number of neighbors is greater than the actual number of neighbors in the edge cases. Thus, this function should only be used to reserve memory for the neighbor_ids vector.
Returns
const DistType The expected number of neighbors

Implements Utopia::Grid< Space >.

1092 {
1093 // Distinguish by mode
1094 if (nb_mode == NBMode::empty) {
1095 return 0;
1096 }
1097 else if (nb_mode == NBMode::Moore) {
1098 // Neighborhood size depends on the distance parameter. Get it.
1099 const auto distance = get_nb_param_distance(nb_params);
1100
1101 if (distance <= 1) {
1102 return std::pow(2 + 1, dim) - 1;
1103 }
1104 else {
1105 return std::pow(2 * distance + 1, dim) - 1;
1106 }
1107 }
1108 else if (nb_mode == NBMode::vonNeumann) {
1109 // Neighborhood size depends on the distance parameter. Get it.
1110 const auto distance = get_nb_param_distance(nb_params);
1111
1112 // This one is more complicated ...
1113 // Define a lambda that can be called recursively
1114 auto num_nbs_impl = [](const unsigned short int d, // dimension
1115 DistType distance,
1116 auto& num_nbs_ref)
1117 {
1118 if (d == 1) {
1119 // End of recursion
1120 return 2 * distance;
1121 }
1122 else {
1123 // Recursive branch
1124 DistType cnt = 0;
1125 while (distance > 0) {
1126 cnt += 2 * num_nbs_ref(d-1, distance, num_nbs_ref);
1127 --distance;
1128 }
1129 return cnt;
1130 }
1131 };
1132
1133 // Call the recursive lambda with the space's dimensionality.
1134 return num_nbs_impl(this->dim, distance, num_nbs_impl);
1135 }
1136 else {
1137 throw std::invalid_argument("No '" + nb_mode_to_string(nb_mode)
1138 + "' available for rectangular grid discretization!");
1139 }
1140 };
const NBMode & nb_mode() const
Const reference to the currently selected neighborhood mode.
Definition base.hh:198
const Config & nb_params() const
The neighborhood parameters of the currently selected neighborhood.
Definition base.hh:203
DistType get_nb_param_distance(const Config &params) const
Extract the distance neighborhood parameter from the given config.
Definition square.hh:1159
std::string nb_mode_to_string(const NBMode &nb_mode)
Given an NBMode enum value, return the corresponding string key.
Definition base.hh:78
@ vonNeumann
The vonNeumann neighborhood, i.e. only nearest neighbors.
@ Moore
The Moore neighborhood, i.e. nearest and next nearest neighbors.
@ empty
Every entity is utterly alone in the world.

◆ extent_of()

template<class Space >
SpaceVec Utopia::SquareGrid< Space >::extent_of ( const IndexType  ) const
inlineoverridevirtual

Returns the extent of the cell with the given ID.

Note
This method does not perform bounds checking of the given ID!

Implements Utopia::Grid< Space >.

149 {
150 return _cell_extent;
151 }

◆ get_nb_func()

template<class Space >
NBFuncID< Base > Utopia::SquareGrid< Space >::get_nb_func ( NBMode  nb_mode,
const Config nb_params 
)
inlineoverrideprotectedvirtual

Retrieve the neighborhood function depending on the mode and parameters.

Implements Utopia::Grid< Space >.

405 {
406 if (nb_mode == NBMode::empty) {
407 return this->_nb_empty;
408 }
409 else if (nb_mode == NBMode::vonNeumann) {
411 }
412 else if (nb_mode == NBMode::Moore) {
414 }
415 else {
416 throw std::invalid_argument("No '" + nb_mode_to_string(nb_mode)
417 + "' available for rectangular grid discretization!");
418 }
419 }
NBFuncID< Self > _nb_empty
A neighborhood function for empty neighborhood.
Definition base.hh:323
NBFuncID< Base > get_nb_func_vonNeumann(const Config &nb_params)
Returns a standalone von-Neumann neighborhood function.
Definition square.hh:441
NBFuncID< Base > get_nb_func_Moore(const Config &nb_params)
Returns a standalone Moore neighborhood function.
Definition square.hh:713

◆ get_nb_func_Moore()

template<class Space >
NBFuncID< Base > Utopia::SquareGrid< Space >::get_nb_func_Moore ( const Config nb_params)
inlineprotected

Returns a standalone Moore neighborhood function.

It extracts the distance parameter from the configuration and depending on the distance parameter and the periodicity of the space decides between four different neighborhood calculation functions. The returned callable does rely on the SquareGrid object, but it includes all parameters from the configuration that can be computed once and then captured; this avoids recomputation.

Parameters
nb_paramsThe configuration for the Moore neighborhood method. Expected keys: distance (optional, defaults to 1), which refers to the Chebychev distance of the included neighbors.
713 {
714 // Extract the optional distance parameter
715 constexpr bool check_shape = true;
716 const auto distance = get_nb_param_distance<check_shape>(nb_params);
717
718 // For distance 1, use the specialized functions which are defined as
719 // class members (to don't bloat this method even more). Those
720 // functions do not require any calculation or capture that goes beyond
721 // the capture of ``this``.
722 if (distance <= 1) {
723 if (this->is_periodic()) {
724 return _nb_Moore_periodic;
725 }
726 else {
728 }
729 }
730 // else: distance is > 1. Depending on periodicity of the grid, define
731 // the relevant lambda and let it capture as many values as possible in
732 // order to avoid recomputation.
733
734 // Compute neighborhood size to allow it to be captured
736 nb_params);
737
738 if (this->is_periodic()) {
739 return [this, distance, nb_size](const IndexType root_id){
740 static_assert(dim == 2,
741 "Moore neighborhood is only available in 2D!");
742
743 // Generate vector in which to store the neighbors...
745
746 // ... and allocate space
747 neighbor_ids.reserve(nb_size);
748
749 // Get all neighbors in the first dimension
750 for (DistType dist=1; dist <= distance; ++dist) {
755 }
756
757 // For these neighbors, add _their_ neighbors in the second dimension
758 for (const auto& nb : neighbor_ids) {
759 for (DistType dist=1; dist <= distance; ++dist) {
764 }
765 }
766
767 // And finally, add the root cell's neighbors in the second dimension
768 for (DistType dist=1; dist <= distance; ++dist) {
773 }
774
775 return neighbor_ids;
776 }; // End of lambda definition for periodic space
777 }
778 else { // Space is non-periodic
779 return [this, distance, nb_size](const IndexType root_id){
780 static_assert(dim == 2,
781 "Moore neighborhood is only available in 2D!");
782
783 // Generate vector in which to store the neighbors...
785
786 // ... and allocate space
787 neighbor_ids.reserve(nb_size);
788
789 // Get all neighbors in the first dimension
790 for (DistType dist=1; dist <= distance; ++dist) {
795 }
796
797 // For these neighbors, add _their_ neighbors in the second dimension
798 for (const auto& nb : neighbor_ids) {
799 for (DistType dist=1; dist <= distance; ++dist) {
803 }
804 }
805
806 // And finally, add the root cell's neighbors in the second dimension
807 for (DistType dist=1; dist <= distance; ++dist) {
812 }
813
814 return neighbor_ids;
815 }; // End of lambda definition for non-periodic space
816 }
817 }
auto nb_size() const
Maximum size of the currently selected neighborhood.
Definition base.hh:208
NBFuncID< Base > _nb_Moore_nonperiodic
Moore neighbors for non-periodic 2D grid.
Definition square.hh:844
NBFuncID< Base > _nb_Moore_periodic
Moore neighbors for periodic 2D grid.
Definition square.hh:820
DistType expected_num_neighbors(const NBMode &nb_mode, const Config &nb_params) const override
Computes the expected number of neighbors for a neighborhood mode.
Definition square.hh:1090
std::vector< IndexType > IndexContainer
Type for container of indices.
Definition types.hh:43
std::size_t IndexType
Type for indices, i.e. values used for container indexing, agent IDs, ...
Definition types.hh:40

◆ get_nb_func_vonNeumann()

template<class Space >
NBFuncID< Base > Utopia::SquareGrid< Space >::get_nb_func_vonNeumann ( const Config nb_params)
inlineprotected

Returns a standalone von-Neumann neighborhood function.

It extracts the distance parameter from the configuration and depending on the distance parameter and the periodicity of the space decides between four different neighborhood calculation functions. The returned callable does rely on the SquareGrid object, but it includes all parameters from the configuration that can be computed once and then captured; this avoids recomputation.

Parameters
nb_paramsThe configuration for the von-Neumann neighborhood method. Expected keys: distance (optional, defaults to 1), which refers to the Manhattan distance of included neighbors.
441 {
442 // Extract the optional distance parameter
443 constexpr bool check_shape = true;
444 const auto distance = get_nb_param_distance<check_shape>(nb_params);
445
446 // For distance 1, use the specialized functions which are defined as
447 // class members (to don't bloat this method even more). Those
448 // functions do not require any calculation or capture that goes beyond
449 // the capture of ``this``.
450 if (distance <= 1) {
451 if (this->is_periodic()) {
453 }
454 else {
456 }
457 }
458 // else: distance is > 1. Depending on periodicity of the grid, define
459 // the relevant lambda and let it capture as many values as possible in
460 // order to avoid recomputation.
461
462 // Compute neighborhood size to allow it to be captured
464 nb_params);
465
466 if (this->is_periodic()) {
467 return [this, distance, nb_size](const IndexType root_id){
468 static_assert((dim >= 1 and dim <= 2),
469 "VonNeumann neighborhood is implemented only for 1D or 2D "
470 "space!");
471
472 // Instantiate container in which to store the neighboring IDs
474
475 // Pre-allocating space brings factor ~2 speed improvement
476 neighbor_ids.reserve(nb_size);
477
478 // Depending on the number of dimensions, add the IDs of
479 // neighboring cells in those dimensions.
480 // Add neighbors in dimension 1
481 for (DistType dist=1; dist <= distance; ++dist) {
483 dist,
486 dist,
488 }
489
490 // If the dimension exists, add neighbors in dimension 2
491 if constexpr (dim >= 2) {
492 // Go through all the previously added neighbors and add
493 // the additional neighbors from the other dimension.
494 // NOTE that this algorithm requires the neighbors nearest
495 // to the root_id to have been pushed to the vector
496 // first. The fixed ordering of the previous addition
497 // is required.
498 const auto nb_size = neighbor_ids.size();
499
500 for (DistType i=0; i < nb_size; ++i) {
501 // Add all neighbor ids up to the maximal distance
502 // along the dimension 2.
503 for (DistType dist = 1;
504 dist <= distance - 1 - i/2 + (i%2)/2;
505 ++dist)
506 {
507 // front neighbor
509 dist,
511
512 // back neighbor
514 dist,
516 }
517 }
518
519 // Finally, add the root cell's neighbors in the 2nd dim.
520 for (DistType dist=1; dist <= distance; ++dist) {
522 dist,
525 dist,
527 }
528 }
529
530 // Return the container of cell indices
531 return neighbor_ids;
532 }; // End of lambda definition for periodic space
533 }
534
535 else { // Space is non-periodic
536 return [this, distance, nb_size](const IndexType root_id){
537 static_assert(((dim == 1) or (dim == 2)),
538 "VonNeumann neighborhood is implemented only for 1D or 2D "
539 "space!");
540
541 // Instantiate containers in which to store the neighboring IDs
544
545 // Pre-allocating space brings a factor ~2 speed improvement
546 // NOTE The front_nb_ids vector needs to reserve memory for all
547 // neighbors including the back neighbors because these
548 // will be added to the container directly before
549 // returning it.
550 front_nb_ids.reserve(nb_size);
551 back_nb_ids.reserve(nb_size / 2);
552
553 // Depending on the number of dimensions, add the IDs of
554 // neighboring cells in those dimensions.
555 // Add front neighbors in dimension 1
556 for (DistType dist=1; dist <= distance; ++dist) {
558 dist,
561 dist,
563 }
564
565 // If the dimension exists, add neighbors in dimension 2
566 if constexpr (dim >= 2) {
567 // Go through the front neighbor ids in dimension 1 and add
568 // the neighbor ids in dimension 2
569 // NOTE that this algorithm requires the neighbors nearest
570 // to the root_id to have been pushed to the vector
571 // first. The fixed ordering of the previous addition
572 // is required.
573 const DistType front_nb_size = front_nb_ids.size();
574
575 for (DistType i=0; i < front_nb_size; ++i) {
576 // Add all front neighbor ids up to the maximal
577 // distance along dimension 2.
578 for (DistType dist = 1;
579 dist <= distance - (i + 1);
580 ++dist)
581 {
582 // add front neighbors ids in dimension 2
584 dist,
586
587 // add back neighbors ids in dimension 2
589 dist,
591 }
592 }
593
594 // Go through the front neighbor ids in dimension 1 and add
595 // the neighbor ids in dimension 2
596 // NOTE that this algorithm requires the neighbors nearest
597 // to the root_id to have been pushed to the vector
598 // first. The fixed ordering of the previous addition
599 // is required.
600 const DistType back_nb_size = back_nb_ids.size();
601
602 for (DistType i=0; i<back_nb_size; ++i) {
603 // Add all neighbor ids up to the maximal distance
604 // along second dimension
605 for (DistType dist = 1;
606 dist <= distance - (i + 1);
607 ++dist)
608 {
610 dist,
612
614 dist,
616 }
617 }
618
619 // Finally, add the root cell's neighbors in the 2nd dim.
620 for (DistType dist=1; dist <= distance; ++dist) {
622 dist,
624
626 dist,
628 }
629 }
630
631 // Combine the front and back neighbors container
632 front_nb_ids.insert(front_nb_ids.end(),
633 back_nb_ids.begin(), back_nb_ids.end());
634
635 // Done now. The front neighbor container contains all the IDs.
636 return front_nb_ids;
637 }; // End of lambda definition for non-periodic space
638 }
639 }
NBFuncID< Base > _nb_vonNeumann_periodic
The Von-Neumann neighborhood for periodic grids.
Definition square.hh:643
NBFuncID< Base > _nb_vonNeumann_nonperiodic
The Von-Neumann neighborhood for non-periodic grids.
Definition square.hh:670

◆ get_nb_param_distance()

template<class Space >
template<bool check_shape = false>
DistType Utopia::SquareGrid< Space >::get_nb_param_distance ( const Config params) const
inlineprotected

Extract the distance neighborhood parameter from the given config.

Note
This will never return a KeyError; if the key is not given, this method will return 1.
Template Parameters
check_shapeWhether to check the shape of the grid is large enough for this distance. For all SquareGrid neighborhoods in periodic space, the grid needs to be at least 2 * distance + 1 cells wide in each dimension.
Parameters
paramsThe neighborhood parameters to extract the distance parameter from.
1159 {
1160 const auto distance = get_as<DistType>("distance", params, 1);
1161
1162 // Check the value is smaller than the grid shape. It needs to fit into
1163 // the shape of the grid, otherwise all the algorithms above would have
1164 // to check for duplicate entries and be set-based, which would be
1165 // very inefficient.
1166 if constexpr (check_shape) {
1167 if ( this->is_periodic()
1168 and (distance * 2 + 1 > this->shape().min()))
1169 {
1170 // To inform about the grid shape, print it to the stringstream
1171 // and include it in the error message below.
1172 std::stringstream shape_ss;
1173 this->shape().print(shape_ss, "Grid Shape:");
1174
1175 throw std::invalid_argument("The grid shape is too small to "
1176 "accomodate a neighborhood with 'distance' parameter set "
1177 "to " + get_as<std::string>("distance", params, "1")
1178 + " in a periodic space!\n" + shape_ss.str());
1179 }
1180 }
1181
1182 return distance;
1183 }

◆ id_shift()

template<class Space >
template<DimType axis>
constexpr IndexType Utopia::SquareGrid< Space >::id_shift ( ) const
inlineconstexprprotected

Return the shift in cell indices necessary if moving along an axis.

It returns in the different cases:

  • axis == 0 -> 1
  • axis == 1 -> shape[0] * 1
  • axis == 2 -> shape[1] * shape[0] * 1
  • axis == 3 -> shape[2] * shape[1] * shape[0] * 1
  • ...
Template Parameters
axisIn which dimension the shift is desired
Returns
constexpr IndexType
891 {
892 if constexpr (axis == 0) {
893 // End of recursion
894 return 1;
895 }
896 else {
897 // Recursive branch
898 return _shape[axis - 1] * id_shift<axis-1>();
899 }
900 }
constexpr IndexType id_shift() const
Return the shift in cell indices necessary if moving along an axis.
Definition square.hh:891

◆ midx_of()

template<class Space >
MultiIndex Utopia::SquareGrid< Space >::midx_of ( const IndexType  id) const
inlineoverridevirtual

Returns the multi-index of the cell with the given ID.

Note
This method does not perform bounds checking of the given ID!

Implements Utopia::Grid< Space >.

125 {
126 static_assert(dim <= 2, "MultiIndex only implemented for 1D and 2D!");
127
128 if constexpr (dim == 1) {
129 return MultiIndex({id % _shape[0]});
130 }
131 else {
132 return MultiIndex({id % _shape[0],
133 id / _shape[0]});
134 }
135 }

◆ num_cells()

template<class Space >
IndexType Utopia::SquareGrid< Space >::num_cells ( ) const
inlineoverridevirtual

Number of square cells required to fill the physical space.

This is calculated simply from the _shape member.

Implements Utopia::Grid< Space >.

96 {
97 return std::accumulate(this->_shape.begin(), this->_shape.end(),
98 1, std::multiplies<IndexType>());
99 };

◆ shape()

template<class Space >
MultiIndex Utopia::SquareGrid< Space >::shape ( ) const
inlineoverridevirtual

Get shape of the square grid.

Implements Utopia::Grid< Space >.

111 {
112 return _shape;
113 }

◆ structure()

template<class Space >
GridStructure Utopia::SquareGrid< Space >::structure ( ) const
inlineoverridevirtual

Structure of the grid.

Implements Utopia::Grid< Space >.

116 {
118 }
@ square
A square lattice grid.

◆ vertices_of()

template<class Space >
std::vector< SpaceVec > Utopia::SquareGrid< Space >::vertices_of ( const IndexType  id) const
inlineoverridevirtual

Returns the vertices of the cell with the given ID.

Only available for 2D currently; the vertices are given in counter-clockwise order, starting with the position of the bottom left-hand vertex of the cell.

Note
This method does not perform bounds checking of the given ID!

Implements Utopia::Grid< Space >.

160 {
161 static_assert(dim == 2,
162 "SquareGrid::vertices_of is only implemented for 2D!");
163
164 std::vector<SpaceVec> vertices{};
165 vertices.reserve(4);
166
167 // NOTE The %-operator performs element-wise multiplication
168 // Counter-clockwise, starting bottom left ...
169 vertices.push_back(midx_of(id) % _cell_extent);
170 vertices.push_back(vertices[0] + _cell_extent % SpaceVec({1., 0.}));
171 vertices.push_back(vertices[0] + _cell_extent);
172 vertices.push_back(vertices[0] + _cell_extent % SpaceVec({0., 1.}));
173
174 return vertices;
175 }
@ vertices
Iterate over vertices.

Member Data Documentation

◆ _cell_extent

template<class Space >
const SpaceVec Utopia::SquareGrid< Space >::_cell_extent
private

The extent of each cell of this square discretization (same for all)

◆ _nb_Moore_nonperiodic

template<class Space >
NBFuncID<Base> Utopia::SquareGrid< Space >::_nb_Moore_nonperiodic
protected

Moore neighbors for non-periodic 2D grid.

846 {
847 static_assert(dim == 2, "Moore neighborhood is only available in 2D!");
848
849 // Generate vector in which to store the neighbors and allocate space
851 neighbor_ids.reserve(8);
852
853 // Get the neighbors in the second dimension
855 // root not at border: have them at indices 0 and 1 now
856 // root at border: less than two neighbors were added
857
858 // Distinguish these two cases
859 if (neighbor_ids.size() == 2) {
860 // Was not at a boundary.
863 }
864 else if (neighbor_ids.size() == 1) {
865 // Was at a front XOR back boundary in dimension 2
867 }
868 // else: was at front AND back boundary (single row of cells in dim 2)
869
870 // Finally, add the root's neighbors in the first dimension
872
873 return neighbor_ids;
874 };

◆ _nb_Moore_periodic

template<class Space >
NBFuncID<Base> Utopia::SquareGrid< Space >::_nb_Moore_periodic
protected
Initial value:

Moore neighbors for periodic 2D grid.

822 {
823 static_assert(dim == 2, "Moore neighborhood is only available in 2D!");
824
825 // Generate vector in which to store the neighbors and allocate space
827 neighbor_ids.reserve(8);
828
829 // Get the neighbors in the second dimension
831 // ...have these neighbors at indices 0 and 1 now.
832
833 // For these neighbors, add _their_ neighbors in the first dimension
836
837 // And finally, add the root cell's neighbors in the first dimension
839
840 return neighbor_ids;
841 };

◆ _nb_vonNeumann_nonperiodic

template<class Space >
NBFuncID<Base> Utopia::SquareGrid< Space >::_nb_vonNeumann_nonperiodic
protected
Initial value:
=
[this](const IndexType root_id)
{
static_assert(((dim == 1) or (dim == 2)),
"VonNeumann neighborhood is only implemented in 1 or 2 dimensions "
"in space!");
neighbor_ids.reserve(2 * dim);
if constexpr (dim >= 2) {
}
return neighbor_ids;
}

The Von-Neumann neighborhood for non-periodic grids.

672 {
673 static_assert(((dim == 1) or (dim == 2)),
674 "VonNeumann neighborhood is only implemented in 1 or 2 dimensions "
675 "in space!");
676
677 // Instantiate container in which to store the neighboring cell IDs
679
680 // The number of neighbors is known; pre-allocating space brings a
681 // speed improvement of about factor 2
682 neighbor_ids.reserve(2 * dim);
683
684 // Depending on the number of dimensions, add the IDs of neighboring
685 // cells in those dimensions
687
688 if constexpr (dim >= 2) {
690 }
691
692 // Return the container of cell indices
693 return neighbor_ids;
694 };

◆ _nb_vonNeumann_periodic

template<class Space >
NBFuncID<Base> Utopia::SquareGrid< Space >::_nb_vonNeumann_periodic
protected
Initial value:
=
[this](const IndexType root_id)
{
static_assert((dim >= 1 and dim <= 2),
"VonNeumann neighborhood is implemented only for 1D or 2D space!");
neighbor_ids.reserve(2 * dim);
if constexpr (dim >= 2) {
}
return neighbor_ids;
}

The Von-Neumann neighborhood for periodic grids.

645 {
646 static_assert((dim >= 1 and dim <= 2),
647 "VonNeumann neighborhood is implemented only for 1D or 2D space!");
648
649 // Instantiate container in which to store the neighboring cell IDs
651
652 // The number of neighbors is known; pre-allocating space brings a
653 // speed improvement of about factor 2.
654 neighbor_ids.reserve(2 * dim);
655
656 // Depending on the number of dimensions, add the IDs of neighboring
657 // cells in those dimensions
659
660 if constexpr (dim >= 2) {
662 }
663
664 // Return the container of cell indices
665 return neighbor_ids;
666 };

◆ _shape

template<class Space >
const MultiIndex Utopia::SquareGrid< Space >::_shape
private

The (multi-index) shape of the grid, resulting from resolution.

◆ dim

template<class Space >
constexpr DimType Utopia::SquareGrid< Space >::dim = Space::dim
staticconstexpr

The dimensionality of the space to be discretized (for easier access)


The documentation for this class was generated from the following file: