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::HexagonalGrid< Space > Class Template Reference

A grid discretization using hexagonal cells. More...

#include <hexagonal.hh>

Inheritance diagram for Utopia::HexagonalGrid< Space >:
Inheritance graph
[legend]
Collaboration diagram for Utopia::HexagonalGrid< 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

 HexagonalGrid (std::shared_ptr< Space > space, const Config &cfg)
 Construct a hexagonal grid discretization.
 
IndexType num_cells () const override
 Number of hexagonal 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 hexagonal 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.
 
DistType expected_num_neighbors (const NBMode &nb_mode, const Config &) const override
 Computes the expected number of neighbors for a neighborhood mode.
 
NBFuncID< Baseget_nb_func_hexagonal (const Config &nb_params)
 Returns a standalone hexagonal neighborhood function.
 
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<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_hexagonal_periodic
 The Von-Neumann neighborhood for periodic grids.
 
NBFuncID< Base_nb_hexagonal_nonperiodic
 The Von-Neumann neighborhood for non-periodic grids.
 
- 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
 Get shape of the hexagonal grid.
 

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 discretization (same for all)
 

Detailed Description

template<class Space>
class Utopia::HexagonalGrid< Space >

A grid discretization using hexagonal cells.

This is a grid discretization using hexagonal cells

The hexagonal tiling is constructed in pointy-topped orientation (one vertex of the hexagon pointing up) with an even row offset (every even row is displaced half a cell width to the right).

The hexagonal cellular arrangement covers a periodic space perfectly. In this setup the right-most cell in offset rows is half cut, covering the remaining space at the left boundary. Similarly, the tip of the cells in the upper-most row covers the free space at the lower boundary. In non-periodic space these cells are cut and at the opposite boundary space remains under-represented.

  1. Top: pointy tops are cut off
  2. Right: offset row cells are cut in half
  3. Bottom: no cutoff, but unrepresented area the same shape as the top-side cutoff. Points within these triangles are mapped to nearest cell
  4. Left: no cutoff, but unrepresented area the same shape as the right-side cutoff. Points within these half-cells are mapped to the nearest cell.
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

For more information on hexagonal grids, see this tutorial.

This class has been implemented following the above tutorial written by Amit Patel and published on www.redblobgames.com (last accessed version from May 2020).

Member Typedef Documentation

◆ Base

Base class type.

◆ Config

The configuration type.

◆ MultiIndex

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

◆ HexagonalGrid()

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

Construct a hexagonal grid discretization.

Parameters
spaceThe space to construct the discretization for
cfgFurther configuration parameters
85 :
86 Base(space, cfg),
89 {
90 // Make sure the cells really are hexagonal
91 static_assert(dim == 2, "Hexagonal grid is only implemented for 2D "
92 "space!");
93
95
96 // check that the cell extent corresponds to that of pointy-topped
97 // hexagonal cells
98 // NOTE the test only requires the aspect ratio is precise to 2% for
99 // easy user experience. The hexagonal grid requires a space
100 // with an aspect ratio involving a sqrt(3)
101 if ( fabs(_cell_extent[0] / _cell_extent[1] - sqrt(3.) / 2.)
102 > get_as<double>("aspect_ratio_tolerance", cfg, 0.02))
103 {
104 SpaceVec required_space = SpaceVec({1., 0.75 * 2. / sqrt(3.)});
105 required_space[0] = this->_space->extent[0];
106 required_space[1] = _shape[1] * _cell_extent[0] * 1.5 /sqrt(3.);
107 std::stringstream required_space_ss;
109
110 throw std::invalid_argument(fmt::format(
111 "Given the extent of the physical space and the specified "
112 "resolution, a mapping with hexagonal cells could not be "
113 "found! Either adjust the the extent of physical space or the "
114 "resolution of the grid. Alternatively increase the tolerance "
115 "to distorted hexagons or choose another grid. \n"
116 "The required aspect ratio of sqrt(3) / 2 is violated by {} "
117 "(> `aspect_ratio_tolerance` = {})! \n"
118 "With the given resolution, set the space extent to : \n"
119 "{} or increase the resolution.",
120 fabs(_cell_extent[0] / _cell_extent[1] - sqrt(3.) / 2.),
121 get_as<double>("aspect_ratio_tolerance", cfg, 0.02),
123 ));
124
125 }
126 }
const std::shared_ptr< Space > _space
The space that is to be discretized.
Definition base.hh:120
const std::shared_ptr< Space > & space() const
Const reference to the space this grid maps to.
Definition base.hh:299
typename Space::SpaceVec SpaceVec
The type of vectors that have a relation to physical space.
Definition hexagonal.hh:59
const SpaceVec _cell_extent
The extent of each cell of this discretization (same for all)
Definition hexagonal.hh:76
MultiIndex determine_shape() const
Get shape of the hexagonal grid.
Definition hexagonal.hh:413
Grid< Space > Base
Base class type.
Definition hexagonal.hh:53
SpaceVec effective_resolution() const override
The effective cell resolution into each physical space dimension.
Definition hexagonal.hh:140
static constexpr DimType dim
The dimensionality of the space to be discretized (for easier access)
Definition hexagonal.hh:56
const MultiIndex _shape
The (multi-index) shape of the grid, resulting from resolution.
Definition hexagonal.hh:71
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_neighbors_in_()

template<class Space >
template<DimType axis, bool periodic>
void Utopia::HexagonalGrid< 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
587 {
588 // Assure the number of dimensions is supported
589 static_assert(axis < 3);
590
591 // left and right
592 if constexpr (axis == 0) {
593 // Compute a normalized id
594 const IndexType nrm_id = root_id % _shape[0];
595
596 // Check if at low value boundary
597 if (nrm_id == 0) {
598 // left most column
599 if constexpr (periodic) {
600 neighbor_ids.push_back(root_id - 1 + _shape[0]);
601 }
602 // else: not periodic; nothing to add here
603 }
604 else {
605 // Not at boundary; no need for the correction term
606 neighbor_ids.push_back(root_id - 1);
607 }
608
609 // Check if at high value boundary
610 if (nrm_id == _shape[0] - 1) {
611 // right most column
612 if constexpr (periodic) {
613 neighbor_ids.push_back(root_id + 1 - _shape[0]);
614 }
615 // else: not periodic; nothing to add here
616 }
617 else {
618 // Not at boundary; no need for the correction term
619 neighbor_ids.push_back(root_id + 1);
620 }
621 }
622 // top-left and bottom-right
623 else if constexpr (axis == 1) {
624 // Compute normalized ids
625 const IndexType nrm_id_0 = root_id % _shape[0]; // column
626 const IndexType nrm_id_1 = ( (root_id % (_shape[0] * _shape[1]))
627 / _shape[0]); // row
628
629 // Check if at low value boundary
630 if (nrm_id_1 == 0) {
631 // bottom row
632 if constexpr (periodic) {
633 if (nrm_id_0 < _shape[0] - 1) {
634 neighbor_ids.push_back( root_id - _shape[0] + 1
635 + _shape[0] * _shape[1]);
636 }
637 else {
638 // also right most
639 neighbor_ids.push_back( root_id - 2*_shape[0] + 1
640 + _shape[0] * _shape[1]);
641 }
642 }
643 // else: not periodic; nothing to add here
644 }
645 else if (nrm_id_0 == _shape[0] - 1 and nrm_id_1 % 2 == 0) {
646 // right column, offset rows
647 if constexpr (periodic) {
648 neighbor_ids.push_back(root_id - 2*_shape[0] + 1);
649 }
650 // else: not periodic; nothing to add here
651 }
652 else if (nrm_id_1 % 2 == 0) {
653 // offset row
654 // Not at boundary; no need for the correction term
655 neighbor_ids.push_back(root_id - _shape[0] + 1);
656 }
657 else {
658 // non-offset row
659 neighbor_ids.push_back(root_id - _shape[0]);
660 }
661
662 // Check if at high value boundary
663 if (nrm_id_1 == _shape[1] - 1) {
664 // top row
665 if constexpr (periodic) {
666 if (nrm_id_0 > 0) {
667 // is an non-offset row
668 neighbor_ids.push_back( root_id + _shape[0] - 1
669 - _shape[0] * _shape[1]);
670 }
671 else {
672 // left most cell
673 neighbor_ids.push_back( root_id + 2*_shape[0] - 1
674 - _shape[0] * _shape[1]);
675 }
676 }
677 // else: not periodic; nothing to add here
678 }
679 else if (nrm_id_0 == 0 and nrm_id_1 % 2 == 1){
680 // left column impair row
681 if constexpr (periodic) {
682 neighbor_ids.push_back(root_id + 2*_shape[0] - 1);
683 }
684 }
685 else if (nrm_id_1 % 2 == 0) {
686 // Not at boundary; no need for the correction term
687 neighbor_ids.push_back(root_id + _shape[0]);
688 }
689 else {
690 // non-offset row
691 neighbor_ids.push_back(root_id + _shape[0] - 1);
692 }
693 }
694 // top right and bottom left
695 else {
696 // Compute normalized ids
697 const IndexType nrm_id_0 = root_id % _shape[0]; // column
698 const IndexType nrm_id_1 = ( (root_id % (_shape[0] * _shape[1]))
699 / _shape[0]); // row
700
701 if (nrm_id_1 == 0) {
702 // bottom row
703 if constexpr (periodic) {
704 neighbor_ids.push_back( root_id - _shape[0]
705 + _shape[0] * _shape[1]);
706 }
707 // else: not periodic; nothing to add here
708 }
709 else if (nrm_id_0 == 0 and nrm_id_1 % 2 == 1) {
710 // left most, non-offset row
711 if constexpr (periodic) {
712 neighbor_ids.push_back(root_id - 1);
713 }
714 // else: not periodic; nothing to add here
715 }
716 else if (nrm_id_1 % 2 == 1) {
717 neighbor_ids.push_back(root_id - _shape[0] - 1);
718 }
719 else {
720 neighbor_ids.push_back(root_id - _shape[0]);
721 }
722
723 if (nrm_id_1 == _shape[1] - 1) {
724 // top row
725 if constexpr (periodic) {
726 neighbor_ids.push_back( root_id + _shape[0]
727 - _shape[0] * _shape[1]);
728 }
729 // else: not periodic; nothing to add here
730 }
731 else if (nrm_id_0 == _shape[0] - 1 and nrm_id_1 % 2 == 0) {
732 // right column in offset row
733 if constexpr (periodic) {
734 neighbor_ids.push_back(root_id + 1);
735 }
736 // else: not periodic; nothing to add here
737 }
738 else if (nrm_id_1 % 2 == 0) {
739 neighbor_ids.push_back(root_id + _shape[0] + 1);
740 }
741 else {
742 neighbor_ids.push_back(root_id + _shape[0]);
743 }
744 }
745 }
std::size_t IndexType
Type for indices, i.e. values used for container indexing, agent IDs, ...
Definition types.hh:40

◆ barycenter_of()

template<class Space >
SpaceVec Utopia::HexagonalGrid< 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 >.

167 {
168 MultiIndex mid = midx_of(id);
169 if (mid[1] % 2 == 0) {
170 // even row
171 return (mid % SpaceVec({1., 0.75}) + SpaceVec({1., 0.5}))
172 % _cell_extent;
173 }
174 else {
175 // odd row
176 return ( (mid % SpaceVec({1, 0.75}) + SpaceVec({0.5, 0.5}))
177 % _cell_extent);
178 }
179 }
MultiIndexType< dim > MultiIndex
The type of multi-index like arrays, e.g. the grid shape.
Definition hexagonal.hh:62
MultiIndex midx_of(const IndexType id) const override
Returns the multi-index of the cell with the given ID.
Definition hexagonal.hh:160

◆ boundary_cells()

template<class Space >
std::set< IndexType > Utopia::HexagonalGrid< 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: 2D: left, right, bottom, top

Implements Utopia::Grid< Space >.

308 {
309 if ( select != "all"
310 and select != "left" and select != "right"
311 and select != "bottom" and select != "top")
312 {
313 throw std::invalid_argument("Invalid value for argument "
314 "`select` in call to method SquareGrid::boundary_cells! "
315 "Available arguments (for currently selected "
316 "dimensionality) are: "
317 "'all', 'left', 'right', 'bottom', 'top'. Given value: '"
318 + select + "'");
319 }
320
321 // The target set all IDs are to be emplaced in
322 std::set<IndexType> bc_ids;
323
324 // For periodic space, this is easy:
325 if (this->is_periodic()) {
326 return {};
327 }
328
329 // NOTE It is important to use the hinting features of std::set
330 // here, which allow to run the following in amortized
331 // constant time instead of logarithmic with set size.
332 // Below, it always makes sense to hint at inserting right
333 // before the end.
334 // Hint for the first element needs to be the beginning
335 auto hint = bc_ids.begin();
336
337 // Bottom boundary (lowest IDs)
338 if (select == "all" or select == "bottom") {
339 // 0, ..., _shape[0] - 1
340 for (DistType id = 0; id < _shape[0]; id++) {
341 bc_ids.emplace_hint(hint, id);
342 hint = bc_ids.end();
343 }
344 }
345
346 // Left boundary
347 if (select == "left") {
348 // First IDs in _shape[1] rows: 0, _shape[0], 2*_shape[0], ...
349 for (DistType row = 0; row < _shape[1]; row++) {
350 bc_ids.emplace_hint(hint, row * _shape[0]);
351 hint = bc_ids.end();
352 }
353 }
354
355 // Right boundary
356 if (select == "right") {
357 // Last IDs in _shape[1] rows
358 const auto offset = _shape[0] - 1;
359
360 for (DistType row = 0; row < _shape[1]; row++) {
361 bc_ids.emplace_hint(hint, offset + row * _shape[0]);
362 hint = bc_ids.end();
363 }
364 }
365
366 // Left AND right (only for 'all' case, allows better hints)
367 if (select == "all") {
368 // First and last IDs in _shape[1] rows
369 const auto offset = _shape[0] - 1;
370
371 for (DistType row = 0; row < _shape[1]; row++) {
372 // Left boundary cell
373 bc_ids.emplace_hint(hint, row * _shape[0]);
374
375 // Right boundary cell (higher than left cell ID)
376 bc_ids.emplace_hint(bc_ids.end(),
377 offset + row * _shape[0]);
378
379 hint = bc_ids.end();
380 }
381 }
382
383 // Top boundary (highest IDs)
384 if (select == "all" or select == "top") {
385 // _shape[0] * (_shape[1]-1), ..., _shape[0] * _shape[1] - 1
386 for (DistType id = _shape[0] * (_shape[1]-1);
387 id < _shape[0] * _shape[1]; id++)
388 {
389 bc_ids.emplace_hint(hint, id);
390 hint = bc_ids.end();
391 }
392 }
393
394 return bc_ids;
395 }
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::HexagonalGrid< 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.

For non-periodic space only: The offset geometry of the hexagonal lattice does not permit to cover a rectangular non-periodic space perfectly. At the boundary, a position in uncovered space is mapped to the closest cell. Details: The hexagonal cellular arrangement is covering the periodic equivalent perfectly. In this setup the right-most cell in offset rows is half cut, covering the space at the left boundary. Similarly, the tip of the cells in the upper-most row covers the free space at the lower boundary. In non-periodic space these cells are cut and at the opposite boundary space remains under-represented. Positions in this space are mapped to the closest cell.

  1. Top: pointy tops are cut off
  2. Right: offset row cells are cut in half
  3. Bottom: no cutoff, but unrepresented area the same shape as the top-side cutoff. Points within these triangles are mapped to nearest cell
  4. Left: no cutoff, but unrepresented area the same shape as the right-side cutoff. Points within these half-cells are mapped to the nearest cell.
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 >.

248 {
249 SpaceVec ridx = pos;
250 // check position
251 if (this->is_periodic()) {
252 ridx = this->_space->map_into_space(pos);
253 }
254 else if (not this->_space->contains(pos)) {
255 throw std::invalid_argument("The given position is outside "
256 "the non-periodic space associated with this grid!");
257 }
258
259 // relative and centered in cell 0
260 ridx = ridx / _cell_extent - SpaceVec({1., 0.5});
261
262 arma::Col<int>::fixed<dim> midx = {
263 static_cast<int>(std::round(ridx[0] + 2. / 3. * ridx[1])),
264 static_cast<int>(std::round(4. / 3. * ridx[1]))
265 };
266
267 // correct for offset
268 midx[0] -= std::floor(midx[1] / 2.);
269
270 if (not this->is_periodic()) {
271 // remap not represented space to first cell in this row
272 // NOTE this distorts the cells, but the rectangular space is
273 // correctly represented
274 if (midx[0] == -1) {
275 midx[0]++;
276 }
277 if (midx[1] == -1) {
278 midx[1]++;
279 }
280
281 // Associate points on high-value boundaries with boundary cells
282 if (midx[0] == static_cast<int>(_shape[0])) {
283 midx[0]--;
284 }
285 }
286 else {
287 midx[0] = (midx[0] + _shape[0]) % _shape[0];
288 midx[1] = (midx[1] + _shape[1]) % _shape[1];
289 }
290
291 // From the multi index, calculate the corresponding ID
292 return static_cast<IndexType>(midx[0] + (midx[1] * _shape[0]));
293 // Equivalent to:
294 // midx[0] * id_shift<0>()
295 // + midx[1] * id_shift<1>()
296 }

◆ determine_shape()

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

Get shape of the hexagonal grid.

Integer rounding takes place here. A physical space of extents of 2.1 length units in each dimension (resp. the shape of the hexagon) 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.

Note
In non-periodic space the the number of rows needs to be a pair value in pointy-top orientation.
Note the offset rows are cut and the right end of the domain (if not periodic). Therefore the width of the cells is such, that N * width = Length, that is non-offset rows span the entire width.
413 {
415
416 // obtain the sidelength of a unit area hexagon
417 // A = 3^1.5 / 2 s^2
418 double s = sqrt(2. / (3 * sqrt(3)));
419 double width = sqrt(3) * s;
420 double height = 2 * s;
421
422 shape[0] = std::round( this->_space->extent[0] / width
423 * this->_resolution);
424 shape[1] = std::round( this->_space->extent[1] / (0.75 * height)
425 * this->_resolution);
426
427 // in non periodic space pair number of rows required
428 if (this->is_periodic()) {
429 shape[1] = shape[1] - shape[1] % 2;
430 }
431
432 return shape;
433 }
MultiIndex shape() const override
Get shape of the hexagonal grid.
Definition hexagonal.hh:146

◆ effective_resolution()

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

The effective cell resolution into each physical space dimension.

Implements Utopia::Grid< Space >.

140 {
141 return SpaceVec({static_cast<double>(_shape[0]), 0.75*_shape[1]})
142 / this->_space->extent;
143 }

◆ expected_num_neighbors()

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

Computes the expected number of neighbors for a neighborhood mode.

Implements Utopia::Grid< Space >.

458 {
459 if (nb_mode == NBMode::empty) {
460 return 0;
461 }
462 else if (nb_mode == NBMode::hexagonal) {
463 return 6;
464 }
465 else {
466 throw std::invalid_argument("No '" + nb_mode_to_string(nb_mode)
467 + "' neighborhood available for HexagonalGrid! "
468 "Available modes: empty, hexagonal.");
469 }
470 }
const NBMode & nb_mode() const
Const reference to the currently selected neighborhood mode.
Definition base.hh:198
std::string nb_mode_to_string(const NBMode &nb_mode)
Given an NBMode enum value, return the corresponding string key.
Definition base.hh:78
@ hexagonal
The hexagonal neighbourhood, i.e. the neighbourhood on a hexagonal grid.
@ empty
Every entity is utterly alone in the world.

◆ extent_of()

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

Returns the extent of the cell with the given ID.

The cell's extent corresponds to the tip-to-tip distance of the cell

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

Implements Utopia::Grid< Space >.

186 {
187 return _cell_extent;
188 }

◆ get_nb_func()

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

Retrieve the neighborhood function depending on the mode.

Implements Utopia::Grid< Space >.

441 {
442 if (nb_mode == NBMode::empty) {
443 return this->_nb_empty;
444 }
445 else if (nb_mode == NBMode::hexagonal) {
447 }
448 else {
449 throw std::invalid_argument("No '" + nb_mode_to_string(nb_mode)
450 + "' neighborhood available for HexagonalGrid! "
451 "Available modes: empty, hexagonal.");
452 }
453 }
NBFuncID< Self > _nb_empty
A neighborhood function for empty neighborhood.
Definition base.hh:323
const Config & nb_params() const
The neighborhood parameters of the currently selected neighborhood.
Definition base.hh:203
NBFuncID< Base > get_nb_func_hexagonal(const Config &nb_params)
Returns a standalone hexagonal neighborhood function.
Definition hexagonal.hh:492

◆ get_nb_func_hexagonal()

template<class Space >
NBFuncID< Base > Utopia::HexagonalGrid< Space >::get_nb_func_hexagonal ( const Config nb_params)
inlineprotected

Returns a standalone hexagonal 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 hexagonal neighborhood method. Expected keys: distance (optional, defaults to 1), which refers to the Manhattan distance of included neighbors.
492 {
493 // Extract the optional distance parameter
494 constexpr bool check_shape = true;
495 const auto distance = get_nb_param_distance<check_shape>(nb_params);
496
497 // For distance 1, use the specialized functions which are defined as
498 // class members (to don't bloat this method even more). Those
499 // functions do not require any calculation or capture that goes beyond
500 // the capture of ``this``.
501 if (distance <= 1) {
502 if (this->is_periodic()) {
504 }
505 else {
507 }
508 }
509 // else: distance is > 1. Depending on periodicity of the grid, define
510 // the relevant lambda and let it capture as many values as possible in
511 // order to avoid recomputation.
512
513 throw std::invalid_argument(fmt::format(
514 "Hexagonal neighborhood is not implemented for a "
515 "distance larger than 1. Requested distance was {}.",
516 distance
517 ));
518 }
NBFuncID< Base > _nb_hexagonal_nonperiodic
The Von-Neumann neighborhood for non-periodic grids.
Definition hexagonal.hh:544
NBFuncID< Base > _nb_hexagonal_periodic
The Von-Neumann neighborhood for periodic grids.
Definition hexagonal.hh:522

◆ get_nb_param_distance()

template<class Space >
template<bool check_shape = false>
DistType Utopia::HexagonalGrid< 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.
764 {
765 const auto distance = get_as<DistType>("distance", params, 1);
766
767 // Check the value is smaller than the grid shape. It needs to fit into
768 // the shape of the grid, otherwise all the algorithms above would have
769 // to check for duplicate entries and be set-based, which would be
770 // very inefficient.
771 if constexpr (check_shape) {
772 if ( this->is_periodic()
773 and (distance * 2 + 1 > this->shape().min()))
774 {
775 // To inform about the grid shape, print it to the stringstream
776 // and include it in the error message below.
777 std::stringstream shape_ss;
778 this->shape().print(shape_ss, "Grid Shape:");
779
780 throw std::invalid_argument("The grid shape is too small to "
781 "accomodate a neighborhood with 'distance' parameter set "
782 "to " + get_as<std::string>("distance", params, "1")
783 + " in a periodic space!\n" + shape_ss.str());
784 }
785 }
786
787 return distance;
788 }

◆ midx_of()

template<class Space >
MultiIndex Utopia::HexagonalGrid< 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 >.

160 {
161 return MultiIndex({id % _shape[0], id / _shape[0]});
162 }

◆ num_cells()

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

Number of hexagonal cells required to fill the physical space.

This is calculated simply from the _shape member.

Implements Utopia::Grid< Space >.

135 {
136 return _shape[0] * _shape[1];
137 }

◆ shape()

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

Get shape of the hexagonal grid.

Implements Utopia::Grid< Space >.

146 {
147 return _shape;
148 }

◆ structure()

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

Structure of the grid.

Implements Utopia::Grid< Space >.

151 {
153 }
@ hexagonal
A hexagonal lattice grid.

◆ vertices_of()

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

Returns the vertices of the cell with the given ID.

The vertices (pointy-topped arrangement) are given in counter-clockwise order, starting with the position of the bottom left-hand vertex (8 o'clock) of the cell.

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

Implements Utopia::Grid< Space >.

197 {
198 std::vector<SpaceVec> vertices{};
199 vertices.reserve(6);
200
201 const SpaceVec center = barycenter_of(id);
202
203 vertices.push_back(center + SpaceVec({-0.5, -0.25}) % _cell_extent);
204 vertices.push_back(center + SpaceVec({ 0. , -0.5 }) % _cell_extent);
205 vertices.push_back(center + SpaceVec({ 0.5, -0.25}) % _cell_extent);
206 vertices.push_back(center + SpaceVec({ 0.5, 0.25}) % _cell_extent);
207 vertices.push_back(center + SpaceVec({ 0. , 0.5 }) % _cell_extent);
208 vertices.push_back(center + SpaceVec({-0.5, 0.25}) % _cell_extent);
209
210 return vertices;
211 }
SpaceVec barycenter_of(const IndexType id) const override
Returns the barycenter of the cell with the given ID.
Definition hexagonal.hh:167
@ vertices
Iterate over vertices.

Member Data Documentation

◆ _cell_extent

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

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

The cell's extent corresponds to the tip-to-tip distance

◆ _nb_hexagonal_nonperiodic

template<class Space >
NBFuncID<Base> Utopia::HexagonalGrid< Space >::_nb_hexagonal_nonperiodic
protected
Initial value:
=
[this](const IndexType root_id)
{
neighbor_ids.reserve(2 * dim);
return neighbor_ids;
}
std::vector< IndexType > IndexContainer
Type for container of indices.
Definition types.hh:43

The Von-Neumann neighborhood for non-periodic grids.

546 {
547 // Instantiate container in which to store the neighboring cell IDs
549
550 // The number of neighbors is known; pre-allocating space brings a
551 // speed improvement of about factor 2
552 neighbor_ids.reserve(2 * dim);
553
554 // Depending on the number of dimensions, add the IDs of neighboring
555 // cells in those dimensions
559
560 // Return the container of cell indices
561 return neighbor_ids;
562 };

◆ _nb_hexagonal_periodic

template<class Space >
NBFuncID<Base> Utopia::HexagonalGrid< Space >::_nb_hexagonal_periodic
protected
Initial value:

The Von-Neumann neighborhood for periodic grids.

524 {
525 // Instantiate container in which to store the neighboring cell IDs
527
528 // The number of neighbors is known; pre-allocating space brings a
529 // speed improvement of about factor 2.
530 neighbor_ids.reserve(3 * dim);
531
532 // Depending on the number of dimensions, add the IDs of neighboring
533 // cells in those dimensions
537
538 // Return the container of cell indices
539 return neighbor_ids;
540 };

◆ _shape

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

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

◆ dim

template<class Space >
constexpr DimType Utopia::HexagonalGrid< 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: