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 | Private Types | Private Member Functions | Static Private Member Functions | Private Attributes | Static Private Attributes | Friends | List of all members
internal::hp::DoFLevel Class Reference

#include <deal.II/hp/dof_level.h>

Public Member Functions

void set_dof_index (const unsigned int obj_index, const unsigned int fe_index, const unsigned int local_index, const types::global_dof_index global_index)
 
types::global_dof_index get_dof_index (const unsigned int obj_index, const unsigned int fe_index, const unsigned int local_index) const
 
unsigned int active_fe_index (const unsigned int obj_index) const
 
bool fe_index_is_active (const unsigned int obj_index, const unsigned int fe_index) const
 
void set_active_fe_index (const unsigned int obj_index, const unsigned int fe_index)
 
unsigned int future_fe_index (const unsigned int obj_index) const
 
void set_future_fe_index (const unsigned int obj_index, const unsigned int fe_index)
 
bool future_fe_index_set (const unsigned int obj_index) const
 
void clear_future_fe_index (const unsigned int obj_index)
 
const types::global_dof_indexget_cell_cache_start (const unsigned int obj_index, const unsigned int dofs_per_cell) const
 
std::size_t memory_consumption () const
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 

Private Types

using offset_type = unsigned int
 
using active_fe_index_type = unsigned short int
 
using signed_active_fe_index_type = signed short int
 

Private Member Functions

template<int dim, int spacedim>
void compress_data (const ::hp::FECollection< dim, spacedim > &fe_collection)
 
template<int dim, int spacedim>
void uncompress_data (const ::hp::FECollection< dim, spacedim > &fe_collection)
 
void normalize_active_fe_indices ()
 

Static Private Member Functions

static bool is_compressed_entry (const active_fe_index_type active_fe_index)
 
static active_fe_index_type get_toggled_compression_state (const active_fe_index_type active_fe_index)
 

Private Attributes

std::vector< active_fe_index_typeactive_fe_indices
 
std::vector< active_fe_index_typefuture_fe_indices
 
std::vector< offset_typedof_offsets
 
std::vector< types::global_dof_indexdof_indices
 
std::vector< offset_typecell_cache_offsets
 
std::vector< types::global_dof_indexcell_dof_indices_cache
 

Static Private Attributes

static const active_fe_index_type invalid_active_fe_index
 

Friends

template<int , int >
class ::hp::DoFHandler
 
struct ::internal::hp::DoFHandlerImplementation::Implementation
 
struct ::internal::DoFCellAccessorImplementation::Implementation
 

Detailed Description

This is the class that stores the degrees of freedom on cells in a hp hierarchy. Compared to faces and edges, the task here is simple since each cell can only have a single active finite element index. Consequently, all we need is one long array with DoF indices and one array of offsets where each cell's indices start within the array of indices. This is in contrast to the DoFObjects class where each face or edge may have more than one associated finite element with corresponding degrees of freedom.

The data stored here is represented by three arrays - The active_fe_indices array stores for each cell which finite element is used on this cell. Since some cells are not active on the current level, some entries in this array may represent an invalid value. - The dof_indices array stores for each active cell on the current level the dofs that associated with the interior of the cell, i.e., the dofs_per_line dofs associated with the line in 1d, and dofs_per_quad and dofs_per_hex in 2d and 3d. These numbers are in general smaller than dofs_per_cell. - The dof_offsets array stores, for each cell, the starting point of the dof indices corresponding to this cell in the dof_indices array. This is analogous to how we store data in compressed row storage for sparse matrices. For cells that are not active on the current level, we store an invalid value for the starting index.

Compression

It is common for the indices stored in dof_indices for one cell to be numbered consecutively. For example, using the standard numbering (without renumbering DoFs), the quad dofs on the first cell of a mesh when using a \(Q_3\) element will be numbered 12, 13, 14, 15. This allows for compression if we only store the first entry and have some way to mark the DoFs on this object as compressed. Here, compression means that we know that subsequent DoF indices can be obtained from the previous ones by just incrementing them by one – in other words, we use a variant of doing run-length encoding. The way to do this is that we use positive FE indices for uncompressed sets of DoFs and if a set of indices is compressed, then we instead store the FE index in binary complement (which we can identify by looking at the sign bit when interpreting the number as a signed one). There are two functions, compress_data() and uncompress_data() that convert between the two possible representations.

Note that compression is not always possible. For example, if one renumbered the example above using DoFRenumbering::downstream with \((1,0)^T\) as direction, then they would likely be numbered 12, 14, 13, 15, which can not be compressed using run-length encoding.

Definition at line 109 of file dof_level.h.

Member Typedef Documentation

◆ offset_type

using internal::hp::DoFLevel::offset_type = unsigned int
private

The type in which we store the offsets into the dof_indices array.

Definition at line 115 of file dof_level.h.

◆ active_fe_index_type

using internal::hp::DoFLevel::active_fe_index_type = unsigned short int
private

The type in which we store the active FE index.

Definition at line 120 of file dof_level.h.

◆ signed_active_fe_index_type

A signed type that matches the type in which we store the active FE index. We use this in computing binary complements.

Definition at line 126 of file dof_level.h.

Member Function Documentation

◆ is_compressed_entry()

bool internal::hp::DoFLevel::is_compressed_entry ( const active_fe_index_type  active_fe_index)
inlinestaticprivate

Given an active_fe_index, return whether the corresponding set of DoF indices are compressed. See the general documentation of this class for a description of when this is the case.

Definition at line 386 of file dof_level.h.

◆ get_toggled_compression_state()

DoFLevel::active_fe_index_type internal::hp::DoFLevel::get_toggled_compression_state ( const active_fe_index_type  active_fe_index)
inlinestaticprivate

Given an active_fe_index (either corresponding to an uncompressed or compressed state), return the active_fe_index that corresponds to the respectively other state. See the general documentation of this class for a description of how compression is indicated.

Definition at line 394 of file dof_level.h.

◆ set_dof_index()

void internal::hp::DoFLevel::set_dof_index ( const unsigned int  obj_index,
const unsigned int  fe_index,
const unsigned int  local_index,
const types::global_dof_index  global_index 
)
inline

Set the global index of the local_index-th degree of freedom located on the object with number obj_index to the value given by global_index. The dof_handler argument is used to access the finite element that is to be used to compute the location where this data is stored.

The third argument, fe_index, denotes which of the finite elements associated with this object we shall access. Refer to the general documentation of the internal::hp::DoFLevel class template for more information.

Definition at line 437 of file dof_level.h.

◆ get_dof_index()

types::global_dof_index internal::hp::DoFLevel::get_dof_index ( const unsigned int  obj_index,
const unsigned int  fe_index,
const unsigned int  local_index 
) const
inline

Return the global index of the local_index-th degree of freedom located on the object with number obj_index. The dof_handler argument is used to access the finite element that is to be used to compute the location where this data is stored.

The third argument, fe_index, denotes which of the finite elements associated with this object we shall access. Refer to the general documentation of the internal::hp::DoFLevel class template for more information.

Definition at line 406 of file dof_level.h.

◆ active_fe_index()

unsigned int internal::hp::DoFLevel::active_fe_index ( const unsigned int  obj_index) const
inline

Return the fe_index of the active finite element on this object.

Definition at line 464 of file dof_level.h.

◆ fe_index_is_active()

bool internal::hp::DoFLevel::fe_index_is_active ( const unsigned int  obj_index,
const unsigned int  fe_index 
) const
inline

Check whether a given finite element index is used on the present object or not.

Definition at line 477 of file dof_level.h.

◆ set_active_fe_index()

void internal::hp::DoFLevel::set_active_fe_index ( const unsigned int  obj_index,
const unsigned int  fe_index 
)
inline

Set the fe_index of the active finite element on this object.

Definition at line 486 of file dof_level.h.

◆ future_fe_index()

unsigned int internal::hp::DoFLevel::future_fe_index ( const unsigned int  obj_index) const
inline

Return the fe_index of the future finite element on this object. If no future_fe_index has been specified, return the active_fe_index instead.

Definition at line 513 of file dof_level.h.

◆ set_future_fe_index()

void internal::hp::DoFLevel::set_future_fe_index ( const unsigned int  obj_index,
const unsigned int  fe_index 
)
inline

Set the fe_index of the future finite element on this object.

Definition at line 526 of file dof_level.h.

◆ future_fe_index_set()

bool internal::hp::DoFLevel::future_fe_index_set ( const unsigned int  obj_index) const
inline

Return whether a future fe index has been set on this object.

Definition at line 553 of file dof_level.h.

◆ clear_future_fe_index()

void internal::hp::DoFLevel::clear_future_fe_index ( const unsigned int  obj_index)
inline

Revoke the future finite element assigned to this object.

Definition at line 563 of file dof_level.h.

◆ get_cell_cache_start()

const types::global_dof_index * internal::hp::DoFLevel::get_cell_cache_start ( const unsigned int  obj_index,
const unsigned int  dofs_per_cell 
) const
inline

Return a pointer to the beginning of the DoF indices cache for a given cell.

Parameters
obj_indexThe number of the cell we are looking at.
dofs_per_cellThe number of DoFs per cell for this cell. This is not used for the hp case but necessary to keep the interface the same as for the non-hp case.
Returns
A pointer to the first DoF index for the current cell. The next dofs_per_cell indices are for the current cell.

Definition at line 573 of file dof_level.h.

◆ memory_consumption()

std::size_t internal::hp::DoFLevel::memory_consumption ( ) const

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

Definition at line 250 of file dof_level.cc.

◆ serialize()

template<class Archive >
void internal::hp::DoFLevel::serialize ( Archive &  ar,
const unsigned int  version 
)
inline

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

Definition at line 595 of file dof_level.h.

◆ compress_data()

template<int dim, int spacedim>
template void internal::hp::DoFLevel::compress_data ( const ::hp::FECollection< dim, spacedim > &  fe_collection)
private

Compress the arrays that store dof indices by using a variant of run- length encoding. See the general documentation of this class for more information.

Parameters
fe_collectionThe object that can tell us how many degrees of freedom each of the finite elements has that we store in this object.

Definition at line 36 of file dof_level.cc.

◆ uncompress_data()

template<int dim, int spacedim>
template void internal::hp::DoFLevel::uncompress_data ( const ::hp::FECollection< dim, spacedim > &  fe_collection)
private

Uncompress the arrays that store dof indices by using a variant of run-length encoding. See the general documentation of this class for more information.

Parameters
fe_collectionThe object that can tell us how many degrees of freedom each of the finite elements has that we store in this object.

Definition at line 167 of file dof_level.cc.

◆ normalize_active_fe_indices()

void internal::hp::DoFLevel::normalize_active_fe_indices ( )
private

Restore the active fe indices stored by the current object to an uncompressed state. At the same time, do not touch any of the other data structures. This leaves the active_fe_indices array out-of-synch from the other member variables, and the function can consequently only be used when the remaining data structures are going to be rebuilt next. This is specifically the case in hp::DoFHandler::distribute_dofs() where we throw away all data in the DoFLevels objects except for the active_fe_indices array. This function therefore simply makes sure that that one array is in uncompressed format so that we can use the information about active fe_indices for all cells that were used in the previous mesh refinement cycle (or the previous time distribute_dofs() was called) without having to care about any of the other data fields.

Definition at line 263 of file dof_level.cc.

Friends And Related Function Documentation

◆ ::hp::DoFHandler

template<int , int >
friend class ::hp::DoFHandler
friend

Definition at line 374 of file dof_level.h.

◆ ::internal::hp::DoFHandlerImplementation::Implementation

friend struct ::internal::hp::DoFHandlerImplementation:: Implementation
friend

Definition at line 376 of file dof_level.h.

◆ ::internal::DoFCellAccessorImplementation::Implementation

friend struct ::internal::DoFCellAccessorImplementation:: Implementation
friend

Definition at line 378 of file dof_level.h.

Member Data Documentation

◆ invalid_active_fe_index

const DoFLevel::active_fe_index_type internal::hp::DoFLevel::invalid_active_fe_index
staticprivate
Initial value:
=

Invalid active_fe_index which will be used as a default value to determine whether a future_fe_index has been set or not.

Definition at line 132 of file dof_level.h.

◆ active_fe_indices

std::vector<active_fe_index_type> internal::hp::DoFLevel::active_fe_indices
private

Indices specifying the finite element of hp::FECollection to use for the different cells on the current level. The vector stores one element per cell since the active_fe_index is unique for cells.

Each active_fe_index will be preset to zero, i.e. the default finite element. However, only active cells are eligible to have a finite element assigned, which will be verified by corresponding accessor functions of the DoFCellAccessor class.

Definition at line 161 of file dof_level.h.

◆ future_fe_indices

std::vector<active_fe_index_type> internal::hp::DoFLevel::future_fe_indices
private

Indices specifying the finite element of hp::FECollection that the cell will be assigned to after the triangulation changed (due to refinement and/or coarsening, or repartitioning). The vector stores one element per cell since the future_fe_index is unique for cells.

Each future_fe_index will be preset to invalid_active_fe_index, indicating that the active_fe_index will remain unchanged. However, only active cells are eligible to have a finite element assigned, which will be verified by corresponding accessor functions of the DoFCellAccessor class.

Definition at line 175 of file dof_level.h.

◆ dof_offsets

std::vector<offset_type> internal::hp::DoFLevel::dof_offsets
private

Store the start index for the degrees of freedom of each object in the dof_indices array. If the cell corresponding to a particular index in this array is not active on this level, then we do not store any DoFs for it. In that case, the offset we store here must be an invalid number and indeed we store (std::vector<types::global_dof_index>::size_type)(-1) for it.

The type we store is then obviously the type the dof_indices array uses for indexing.

Definition at line 189 of file dof_level.h.

◆ dof_indices

std::vector<types::global_dof_index> internal::hp::DoFLevel::dof_indices
private

Store the global indices of the degrees of freedom. information. The dof_offsets field determines where each (active) cell's data is stored.

Definition at line 196 of file dof_level.h.

◆ cell_cache_offsets

std::vector<offset_type> internal::hp::DoFLevel::cell_cache_offsets
private

The offsets for each cell of the cache that holds all DoF indices.

Definition at line 201 of file dof_level.h.

◆ cell_dof_indices_cache

std::vector<types::global_dof_index> internal::hp::DoFLevel::cell_dof_indices_cache
private

Cache for the DoF indices on cells. The size of this array equals the sum over all cells of selected_fe[active_fe_index[cell]].dofs_per_cell.

Definition at line 208 of file dof_level.h.


The documentation for this class was generated from the following files:
internal::hp::DoFLevel::active_fe_index_type
unsigned short int active_fe_index_type
Definition: dof_level.h:120