Reference documentation for deal.II version 9.2.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Public Member Functions | Static Public Member Functions | Public Attributes | List of all members
internal::TriangulationImplementation::TriaLevel< dim > Class Template Reference

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

Public Member Functions

void reserve_space (const unsigned int total_cells, const unsigned int dimension, const unsigned int space_dimension)
 
void monitor_memory (const unsigned int true_dimension) const
 
std::size_t memory_consumption () const
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 

Static Public Member Functions

static ::ExceptionBaseExcMemoryInexact (int arg1, int arg2)
 

Public Attributes

std::vector< std::uint8_t > refine_flags
 
std::vector< boolcoarsen_flags
 
std::vector< unsigned intactive_cell_indices
 
std::vector< std::pair< int, int > > neighbors
 
std::vector< types::subdomain_idsubdomain_ids
 
std::vector< types::subdomain_idlevel_subdomain_ids
 
std::vector< intparents
 
std::vector< booldirection_flags
 
TriaObjects< TriaObject< dim > > cells
 

Detailed Description

template<int dim>
class internal::TriangulationImplementation::TriaLevel< dim >

Store all information which belongs to one level of the multilevel hierarchy.

In TriaLevel, all cell data is stored which is not dependent on the dimension, e.g. a field to store the refinement flag for the cells (what a cell actually is is declared elsewhere), etc. See also TriaObjects for non level-oriented data.

There is another field, which may fit in here, namely the material data (for cells) or the boundary indicators (for faces), but since we need for a line or quad either boundary information or material data, we store them with the lines and quads rather than with the common data. Likewise, in 3d, we need boundary indicators for lines and quads (we need to know how to refine a line if the two adjacent faces have different boundary indicators), and material data for cells.

Author
Wolfgang Bangerth, Guido Kanschat, 1998, 2007

Definition at line 58 of file tria_levels.h.

Member Function Documentation

◆ reserve_space()

template<int dim>
void internal::TriangulationImplementation::TriaLevel< dim >::reserve_space ( const unsigned int  total_cells,
const unsigned int  dimension,
const unsigned int  space_dimension 
)

Reserve enough space to accommodate total_cells cells on this level. Since there are no used flags on this level, you have to give the total number of cells, not only the number of newly to accommodate ones, like in the TriaLevel<N>::reserve_space functions, with N>0.

Since the number of neighbors per cell depends on the dimensions, you have to pass that additionally.

Definition at line 29 of file tria_levels.cc.

◆ monitor_memory()

template<int dim>
void internal::TriangulationImplementation::TriaLevel< dim >::monitor_memory ( const unsigned int  true_dimension) const

Check the memory consistency of the different containers. Should only be called with the preprocessor flag DEBUG set. The function should be called from the functions of the higher TriaLevel classes.

Definition at line 91 of file tria_levels.cc.

◆ memory_consumption()

template<int dim>
std::size_t internal::TriangulationImplementation::TriaLevel< dim >::memory_consumption

Determine an estimate for the memory consumption (in bytes) of this object.

Definition at line 103 of file tria_levels.cc.

◆ serialize()

template<int dim>
template<class Archive >
void internal::TriangulationImplementation::TriaLevel< dim >::serialize ( Archive &  ar,
const unsigned int  version 
)

Read or write the data of this object to or from a stream for the purpose of serialization

Definition at line 276 of file tria_levels.h.

Member Data Documentation

◆ refine_flags

template<int dim>
std::vector<std::uint8_t> internal::TriangulationImplementation::TriaLevel< dim >::refine_flags

RefinementCase<dim>::Type flags for the cells to be refined with or not (RefinementCase<dim>::no_refinement). The meaning what a cell is, is dimension specific, therefore also the length of this vector depends on the dimension: in one dimension, the length of this vector equals the length of the lines vector, in two dimensions that of the quads vector, etc.

Definition at line 69 of file tria_levels.h.

◆ coarsen_flags

template<int dim>
std::vector<bool> internal::TriangulationImplementation::TriaLevel< dim >::coarsen_flags

Same meaning as the one above, but specifies whether a cell must be coarsened.

Definition at line 75 of file tria_levels.h.

◆ active_cell_indices

template<int dim>
std::vector<unsigned int> internal::TriangulationImplementation::TriaLevel< dim >::active_cell_indices

An integer that, for every active cell, stores the how many-th active cell this is. For non-active cells, this value is unused and set to an invalid value.

Definition at line 83 of file tria_levels.h.

◆ neighbors

template<int dim>
std::vector<std::pair<int, int> > internal::TriangulationImplementation::TriaLevel< dim >::neighbors

Levels and indices of the neighbors of the cells. Convention is, that the neighbors of the cell with index i are stored in the fields following \(i*(2*real\_space\_dimension)\), e.g. in one spatial dimension, the neighbors of cell 0 are stored in neighbors[0] and neighbors[1], the neighbors of cell 1 are stored in neighbors[2] and neighbors[3], and so on.

In neighbors, neighbors[i].first is the level, while neighbors[i].first is the index of the neighbor.

If a neighbor does not exist (cell is at the boundary), level=index=-1 is set.

Conventions: The ith neighbor of a cell is the one which shares the ith face (Line in 2D, Quad in 3D) of this cell.

The neighbor of a cell has at most the same level as this cell, i.e. it may or may not be refined.

In one dimension, a neighbor may have any level less or equal the level of this cell. If it has the same level, it may be refined an arbitrary number of times, but the neighbor pointer still points to the cell on the same level, while the neighbors of the children of the neighbor may point to this cell or its children.

In two and more dimensions, the neighbor is either on the same level and refined (in which case its children have neighbor pointers to this cell or its direct children), unrefined on the same level or one level down (in which case its neighbor pointer points to the mother cell of this cell).

Definition at line 118 of file tria_levels.h.

◆ subdomain_ids

template<int dim>
std::vector<types::subdomain_id> internal::TriangulationImplementation::TriaLevel< dim >::subdomain_ids

One integer per cell to store which subdomain it belongs to. This field is most often used in parallel computations, where it denotes which processor shall work on/owns the cells with a given subdomain number.

This number is only used on active cells.

Definition at line 128 of file tria_levels.h.

◆ level_subdomain_ids

template<int dim>
std::vector<types::subdomain_id> internal::TriangulationImplementation::TriaLevel< dim >::level_subdomain_ids

The subdomain id used on each level for parallel multigrid.

In contrast to the subdomain_id, this number is also used on inactive cells once the mesh has been partitioned also on the lower levels of the multigrid hierarchy.

Definition at line 137 of file tria_levels.h.

◆ parents

template<int dim>
std::vector<int> internal::TriangulationImplementation::TriaLevel< dim >::parents

One integer for every consecutive pair of cells to store which index their parent has.

(We store this information once for each pair of cells since every refinement, isotropic or anisotropic, and in any space dimension, always creates children in multiples of two, so there is no need to store the parent index for every cell.)

Definition at line 148 of file tria_levels.h.

◆ direction_flags

template<int dim>
std::vector<bool> internal::TriangulationImplementation::TriaLevel< dim >::direction_flags

One bool per cell to indicate the direction of the normal true: use orientation from vertex false: revert the orientation. See GlossDirectionFlag.

This is only used for codim==1 meshes.

Definition at line 157 of file tria_levels.h.

◆ cells

The object containing the data on lines and related functions

Definition at line 162 of file tria_levels.h.


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