Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
Public Types | Public Member Functions | Private Attributes | Friends | List of all members
CellId Class Reference

#include <deal.II/grid/cell_id.h>

Public Types

using binary_type = std::array< unsigned int, 4 >
 

Public Member Functions

 CellId (const unsigned int coarse_cell_id, const std::vector< std::uint8_t > &child_indices)
 
 CellId (const unsigned int coarse_cell_id, const unsigned int n_child_indices, const std::uint8_t *child_indices)
 
 CellId (const binary_type &binary_representation)
 
 CellId ()
 
std::string to_string () const
 
template<int dim>
binary_type to_binary () const
 
template<int dim, int spacedim>
Triangulation< dim, spacedim >::cell_iterator to_cell (const Triangulation< dim, spacedim > &tria) const
 
bool operator== (const CellId &other) const
 
bool operator!= (const CellId &other) const
 
bool operator< (const CellId &other) const
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 

Private Attributes

unsigned int coarse_cell_id
 
unsigned int n_child_indices
 
std::array< std::uint8_t, internal::p4est::functions< 2 >::max_level > child_indices
 

Friends

std::istream & operator>> (std::istream &is, CellId &cid)
 
std::ostream & operator<< (std::ostream &os, const CellId &cid)
 

Detailed Description

A class to represent a unique ID for a cell in a Triangulation. It is returned by cell->id() if cell is a cell iterator.

This class stores the index of the coarse cell from which a cell is descendant, together with information on how to reach the cell from that coarse cell (i.e., which child index to take on each level of the triangulation when moving from one cell to its children). The important point about this class is that an object of the current class uniquely identifies a cell in triangulation, and it even does so in the context of objects of type parallel::distributed::Triangulation where the local portion of a mesh may not store all cells. For example, the CellId computed for a ghost cell on one processor will be exactly the same as the CellId computed for the very same cell on the processor that actually owns the cell, although the level and index of the iterators pointing to that cell within the triangulation stored on each of the processors may (and in general will) be different. In other words, CellId provides the tool with which it is possible to globally and uniquely identify cells in a parallel triangulation, and consequently makes it possible to exchange, between processors, data tied to individual cells.

Note
How this data is internally represented is not of importance (and not exposed on purpose).

Definition at line 63 of file cell_id.h.

Member Typedef Documentation

◆ binary_type

using CellId::binary_type = std::array<unsigned int, 4>

A type that is used to encode the CellId data in a compact and fast way (e.g. for MPI transfer to other processes). Note that it limits the number of children that can be transferred to 20 in 3D and 30 in 2D (using 2 times 32 bit for storage), a limitation that is identical to the one used by p4est.

Definition at line 73 of file cell_id.h.

Constructor & Destructor Documentation

◆ CellId() [1/4]

CellId::CellId ( const unsigned int  coarse_cell_id,
const std::vector< std::uint8_t > &  child_indices 
)

Construct a CellId object with a given coarse_cell_id and vector of child indices. child_indices is interpreted identical to the member variable with the same name, namely each entry denotes which child to pick from one refinement level to the next, starting with the coarse cell, until we get to the cell represented by the current object. Therefore, each entry should be a number between 0 and the number of children of a cell in the current space dimension (i.e., GeometryInfo<dim>::max_children_per_cell).

Definition at line 38 of file cell_id.cc.

◆ CellId() [2/4]

CellId::CellId ( const unsigned int  coarse_cell_id,
const unsigned int  n_child_indices,
const std::uint8_t *  child_indices 
)

Construct a CellId object with a given coarse_cell_id and array of child indices provided in child_indices. child_indices is interpreted identical to the member variable with the same name, namely each entry denotes which child to pick from one refinement level to the next, starting with the coarse cell, until we get to the cell represented by the current object. Therefore, each entry should be a number between 0 and the number of children of a cell in the current space dimension (i.e., GeometryInfo<dim>::max_children_per_cell). The array child_indices must have at least n_child_indices valid entries.

Definition at line 49 of file cell_id.cc.

◆ CellId() [3/4]

CellId::CellId ( const binary_type binary_representation)

Construct a CellId object with a given binary representation that was previously constructed by CellId::to_binary.

Definition at line 61 of file cell_id.cc.

◆ CellId() [4/4]

CellId::CellId ( )

Construct an invalid CellId.

Definition at line 24 of file cell_id.cc.

Member Function Documentation

◆ to_string()

std::string CellId::to_string ( ) const

Return a human readable string representation of this CellId.

Definition at line 148 of file cell_id.cc.

◆ to_binary()

template<int dim>
CellId::binary_type CellId::to_binary ( ) const

Return a compact and fast binary representation of this CellId.

Definition at line 102 of file cell_id.cc.

◆ to_cell()

template<int dim, int spacedim>
Triangulation< dim, spacedim >::cell_iterator CellId::to_cell ( const Triangulation< dim, spacedim > &  tria) const

Return a cell_iterator to the cell represented by this CellId.

Definition at line 159 of file cell_id.cc.

◆ operator==()

bool CellId::operator== ( const CellId other) const
inline

Compare two CellId objects for equality.

Definition at line 261 of file cell_id.h.

◆ operator!=()

bool CellId::operator!= ( const CellId other) const
inline

Compare two CellIds for inequality.

Definition at line 278 of file cell_id.h.

◆ operator<()

bool CellId::operator< ( const CellId other) const
inline

Compare two CellIds with regard to an ordering. The details of this ordering are unspecified except that the operation provides a total ordering among all cells.

Definition at line 286 of file cell_id.h.

◆ serialize()

template<class Archive >
void CellId::serialize ( Archive &  ar,
const unsigned int  version 
)

Boost serialization function

Serialization function

Definition at line 221 of file cell_id.h.

Friends And Related Function Documentation

◆ operator>>

std::istream& operator>> ( std::istream &  is,
CellId cid 
)
friend

Read a CellId object from a stream.

Definition at line 232 of file cell_id.h.

◆ operator<<

std::ostream& operator<< ( std::ostream &  os,
const CellId cid 
)
friend

Write a CellId object into a stream.

Definition at line 202 of file cell_id.h.

Member Data Documentation

◆ coarse_cell_id

unsigned int CellId::coarse_cell_id
private

The number of the coarse cell within whose tree the cell represented by the current object is located.

Definition at line 166 of file cell_id.h.

◆ n_child_indices

unsigned int CellId::n_child_indices
private

The number of child indices stored in the child_indices array. This is equivalent to (level-1) of the current cell.

Definition at line 172 of file cell_id.h.

◆ child_indices

std::array<std::uint8_t, internal::p4est::functions<2>::max_level> CellId::child_indices
private

An array of integers that denotes which child to pick from one refinement level to the next, starting with the coarse cell, until we get to the cell represented by the current object. Only the first n_child_indices entries are used, but we use a statically allocated array instead of a vector of size n_child_indices to speed up creation of this object. If the given dimensions ever become a limitation the array can be extended.

Definition at line 185 of file cell_id.h.


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