Reference documentation for deal.II version 9.0.0
|
#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) |
const types::global_dof_index * | get_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 | |
typedef unsigned int | offset_type |
typedef unsigned short int | active_fe_index_type |
typedef signed short int | signed_active_fe_index_type |
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_type > | active_fe_indices |
std::vector< offset_type > | dof_offsets |
std::vector< types::global_dof_index > | dof_indices |
std::vector< offset_type > | cell_cache_offsets |
std::vector< types::global_dof_index > | cell_dof_indices_cache |
Friends | |
template<int , int > | |
class | ::hp::DoFHandler |
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.
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 103 of file dof_level.h.
|
private |
The type in which we store the offsets into the dof_indices array.
Definition at line 109 of file dof_level.h.
|
private |
The type in which we store the active FE index.
Definition at line 114 of file dof_level.h.
|
private |
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 120 of file dof_level.h.
|
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 330 of file dof_level.h.
|
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 339 of file dof_level.h.
|
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 384 of file dof_level.h.
|
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 351 of file dof_level.h.
|
inline |
Return the fe_index of the active finite element on this object.
Definition at line 412 of file dof_level.h.
|
inline |
Check whether a given finite element index is used on the present object or not.
Definition at line 428 of file dof_level.h.
|
inline |
Set the fe_index of the active finite element on this object.
Definition at line 439 of file dof_level.h.
|
inline |
Return a pointer to the beginning of the DoF indices cache for a given cell.
obj_index | The number of the cell we are looking at. |
dofs_per_cell | The 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. |
Definition at line 462 of file dof_level.h.
std::size_t internal::hp::DoFLevel::memory_consumption | ( | ) | const |
Determine an estimate for the memory consumption (in bytes) of this object.
Definition at line 226 of file dof_level.cc.
|
inline |
Read or write the data of this object to or from a stream for the purpose of serialization
Definition at line 485 of file dof_level.h.
|
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.
fe_collection | The 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 29 of file dof_level.cc.
|
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.
fe_collection | The 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 150 of file dof_level.cc.
|
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 237 of file dof_level.cc.
|
friend |
Make hp::DoFHandler and its auxiliary class a friend since it is the class that needs to create these data structures.
Definition at line 320 of file dof_level.h.
|
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.
If a cell is not active on the level corresponding to the current object (i.e., it has children on higher levels) then it does not have an associated fe index and we store an invalid fe index marker instead.
Definition at line 151 of file dof_level.h.
|
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 165 of file dof_level.h.
|
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 172 of file dof_level.h.
|
private |
The offsets for each cell of the cache that holds all DoF indices.
Definition at line 177 of file dof_level.h.
|
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 184 of file dof_level.h.