Utopia  2
Framework for studying models of complex & adaptive systems.
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. More...
 
using SpaceVec = typename Space::SpaceVec
 The type of vectors that have a relation to physical space. More...
 
using MultiIndex = MultiIndexType< dim >
 The type of multi-index like arrays, e.g. the grid shape. More...
 
using Config = DataIO::Config
 The configuration type. More...
 
- Public Types inherited from Utopia::Grid< Space >
using Self = Grid< Space >
 Type of this class, i.e. the base grid class. More...
 
using SpaceVec = typename Space::SpaceVec
 The type of vectors that have a relation to physical space. More...
 
using MultiIndex = MultiIndexType< dim >
 The type of multi-index like arrays, e.g. the grid shape. More...
 
using Config = DataIO::Config
 The configuration type. More...
 

Public Member Functions

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

Static Public Attributes

static constexpr DimType dim = Space::dim
 The dimensionality of the space to be discretized (for easier access) More...
 
- 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) More...
 

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. More...
 
NBFuncID< Baseget_nb_func_vonNeumann (const Config &nb_params)
 Returns a standalone von-Neumann neighborhood function. More...
 
NBFuncID< Baseget_nb_func_Moore (const Config &nb_params)
 Returns a standalone Moore neighborhood function. More...
 
template<DimType axis>
constexpr IndexType id_shift () const
 Return the shift in cell indices necessary if moving along an axis. More...
 
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. More...
 
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. More...
 
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. More...
 
DistType expected_num_neighbors (const NBMode &nb_mode, const Config &nb_params) const override
 Computes the expected number of neighbors for a neighborhood mode. More...
 
template<bool check_shape = false>
DistType get_nb_param_distance (const Config &params) const
 Extract the distance neighborhood parameter from the given config. More...
 

Protected Attributes

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

Private Member Functions

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

Private Attributes

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

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

template<class Space >
using Utopia::SquareGrid< Space >::Base = Grid<Space>

Base class type.

◆ Config

template<class Space >
using Utopia::SquareGrid< Space >::Config = DataIO::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

template<class Space >
using Utopia::SquareGrid< Space >::SpaceVec = typename Space::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

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

◆ 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

◆ 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

◆ 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 >.

◆ 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 >.

◆ 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 >.

◆ 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.

◆ 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 >.

◆ 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 >.

◆ 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 >.

◆ 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 >.

◆ 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.

◆ 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.

◆ 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.

◆ 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

◆ 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 >.

◆ 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 >.

◆ shape()

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

Get shape of the square grid.

Implements Utopia::Grid< Space >.

◆ structure()

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

Structure of the grid.

Implements Utopia::Grid< Space >.

◆ 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 >.

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.

◆ _nb_Moore_periodic

template<class Space >
NBFuncID<Base> Utopia::SquareGrid< Space >::_nb_Moore_periodic
protected
Initial value:
=
[this](const IndexType root_id)
{
static_assert(dim == 2, "Moore neighborhood is only available in 2D!");
IndexContainer neighbor_ids{};
neighbor_ids.reserve(8);
add_neighbors_in_<1, true>(root_id, neighbor_ids);
add_neighbors_in_<0, true>(neighbor_ids[0], neighbor_ids);
add_neighbors_in_<0, true>(neighbor_ids[1], neighbor_ids);
add_neighbors_in_<0, true>(root_id, neighbor_ids);
return neighbor_ids;
}
static constexpr DimType dim
The dimensionality of the space to be discretized (for easier access)
Definition: square.hh:37
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

Moore neighbors for periodic 2D grid.

◆ _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!");
IndexContainer neighbor_ids{};
neighbor_ids.reserve(2 * dim);
add_neighbors_in_<0, false>(root_id, neighbor_ids);
if constexpr (dim >= 2) {
add_neighbors_in_<1, false>(root_id, neighbor_ids);
}
return neighbor_ids;
}

The Von-Neumann neighborhood for non-periodic grids.

◆ _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!");
IndexContainer neighbor_ids{};
neighbor_ids.reserve(2 * dim);
add_neighbors_in_<0, true>(root_id, neighbor_ids);
if constexpr (dim >= 2) {
add_neighbors_in_<1, true>(root_id, neighbor_ids);
}
return neighbor_ids;
}

The Von-Neumann neighborhood for periodic grids.

◆ _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: