deal.II version 9.7.1
\(\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\}}\)
Loading...
Searching...
No Matches
dof_handler.h File Reference
#include <deal.II/base/config.h>
#include <deal.II/base/exceptions.h>
#include <deal.II/base/function.h>
#include <deal.II/base/index_set.h>
#include <deal.II/base/iterator_range.h>
#include <deal.II/base/observer_pointer.h>
#include <deal.II/base/types.h>
#include <deal.II/dofs/block_info.h>
#include <deal.II/dofs/dof_accessor.h>
#include <deal.II/dofs/dof_faces.h>
#include <deal.II/dofs/dof_iterator_selector.h>
#include <deal.II/dofs/dof_levels.h>
#include <deal.II/dofs/number_cache.h>
#include <deal.II/hp/fe_collection.h>
#include <boost/serialization/split_member.hpp>
#include <boost/signals2/connection.hpp>
#include <map>
#include <memory>
#include <set>
#include <vector>

Go to the source code of this file.

Classes

class  MGVertexDoFs
struct  ActiveFEIndexTransfer

Namespaces

namespace  internal
namespace  internal::hp
namespace  internal::hp::DoFHandlerImplementation

Typedefs

using cell_accessor = typename ActiveSelector::CellAccessor
using face_accessor = typename ActiveSelector::FaceAccessor
using line_iterator = typename ActiveSelector::line_iterator
using active_line_iterator = typename ActiveSelector::active_line_iterator
using quad_iterator = typename ActiveSelector::quad_iterator
using active_quad_iterator = typename ActiveSelector::active_quad_iterator
using hex_iterator = typename ActiveSelector::hex_iterator
using active_hex_iterator = typename ActiveSelector::active_hex_iterator
using active_cell_iterator = typename ActiveSelector::active_cell_iterator
using cell_iterator = typename ActiveSelector::cell_iterator
using face_iterator = typename ActiveSelector::face_iterator
using active_face_iterator = typename ActiveSelector::active_face_iterator
using level_cell_accessor = typename LevelSelector::CellAccessor
using level_face_accessor = typename LevelSelector::FaceAccessor
using level_cell_iterator = typename LevelSelector::cell_iterator
using level_face_iterator = typename LevelSelector::face_iterator
using offset_type = unsigned int

Functions

 DoFHandler ()
 DoFHandler (const Triangulation< dim, spacedim > &tria)
 DoFHandler (const DoFHandler &)=delete
virtual ~DoFHandler () override
DoFHandleroperator= (const DoFHandler &)=delete
void set_active_fe_indices (const std::vector< types::fe_index > &active_fe_indices)
void set_active_fe_indices (const std::vector< unsigned int > &active_fe_indices)
std::vector< types::fe_indexget_active_fe_indices () const
void get_active_fe_indices (std::vector< unsigned int > &active_fe_indices) const
void set_future_fe_indices (const std::vector< types::fe_index > &future_fe_indices)
std::vector< types::fe_indexget_future_fe_indices () const
void reinit (const Triangulation< dim, spacedim > &tria)
void distribute_dofs (const FiniteElement< dim, spacedim > &fe)
void distribute_dofs (const hp::FECollection< dim, spacedim > &fe)
void distribute_mg_dofs ()
bool has_hp_capabilities () const
bool has_level_dofs () const
bool has_active_dofs () const
void initialize_local_block_info ()
void clear ()
void renumber_dofs (const std::vector< types::global_dof_index > &new_numbers)
void renumber_dofs (const unsigned int level, const std::vector< types::global_dof_index > &new_numbers)
unsigned int max_couplings_between_dofs () const
unsigned int max_couplings_between_boundary_dofs () const
types::global_dof_index n_dofs () const
types::global_dof_index n_dofs (const unsigned int level) const
types::global_dof_index n_boundary_dofs () const
template<typename number>
types::global_dof_index n_boundary_dofs (const std::map< types::boundary_id, const Function< spacedim, number > * > &boundary_ids) const
types::global_dof_index n_boundary_dofs (const std::set< types::boundary_id > &boundary_ids) const
const BlockInfoblock_info () const
types::global_dof_index n_locally_owned_dofs () const
const IndexSetlocally_owned_dofs () const
const IndexSetlocally_owned_mg_dofs (const unsigned int level) const
const FiniteElement< dim, spacedim > & get_fe (const types::fe_index index=0) const
const hp::FECollection< dim, spacedim > & get_fe_collection () const
const Triangulation< dim, spacedim > & get_triangulation () const
MPI_Comm get_mpi_communicator () const
MPI_Comm get_communicator () const
void prepare_for_serialization_of_active_fe_indices ()
void deserialize_active_fe_indices ()
virtual std::size_t memory_consumption () const
template<class Archive>
void save (Archive &ar, const unsigned int version) const
template<class Archive>
void load (Archive &ar, const unsigned int version)
template<class Archive>
void serialize (Archive &archive, const unsigned int version)
static ::ExceptionBaseExcNoFESelected ()
static ::ExceptionBaseExcInvalidBoundaryIndicator ()
static ::ExceptionBaseExcInvalidLevel (int arg1)
static ::ExceptionBaseExcNewNumbersNotConsecutive (types::global_dof_index arg1)
static ::ExceptionBaseExcInvalidFEIndex (int arg1, int arg2)
static ::ExceptionBaseExcOnlyAvailableWithHP ()
static ::ExceptionBaseExcNotImplementedWithHP ()
void clear_space ()
void clear_mg_space ()
void setup_policy ()
void connect_to_triangulation_signals ()
void create_active_fe_table ()
void update_active_fe_table ()
void pre_transfer_action ()
void post_transfer_action ()
void pre_distributed_transfer_action ()
void post_distributed_transfer_action ()
template<int dim, int spacedim>
void internal::hp::DoFHandlerImplementation::communicate_future_fe_indices (DoFHandler< dim, spacedim > &dof_handler)
template<int dim, int spacedim = dim>
unsigned int internal::hp::DoFHandlerImplementation::dominated_future_fe_on_children (const typename DoFHandler< dim, spacedim >::cell_iterator &parent)
static ::ExceptionBaseinternal::hp::DoFHandlerImplementation::ExcNoDominatedFiniteElementOnChildren ()
Cell iterator functions
cell_iterator begin (const unsigned int level=0) const
active_cell_iterator begin_active (const unsigned int level=0) const
cell_iterator end () const
cell_iterator end (const unsigned int level) const
active_cell_iterator end_active (const unsigned int level) const
level_cell_iterator begin_mg (const unsigned int level=0) const
level_cell_iterator end_mg (const unsigned int level) const
level_cell_iterator end_mg () const
Cell iterator functions returning ranges of iterators
IteratorRange< cell_iteratorcell_iterators () const
IteratorRange< active_cell_iteratoractive_cell_iterators () const
IteratorRange< level_cell_iteratormg_cell_iterators () const
IteratorRange< cell_iteratorcell_iterators_on_level (const unsigned int level) const
IteratorRange< active_cell_iteratoractive_cell_iterators_on_level (const unsigned int level) const
IteratorRange< level_cell_iteratormg_cell_iterators_on_level (const unsigned int level) const

Variables

static constexpr unsigned int dimension = dim
static constexpr unsigned int space_dimension = spacedim
static const types::fe_index default_fe_index = 0
BlockInfo block_info_object
bool hp_capability_enabled
ObserverPointer< const Triangulation< dim, spacedim >, DoFHandler< dim, spacedim > > tria
hp::FECollection< dim, spacedim > fe_collection
std::unique_ptr<::internal::DoFHandlerImplementation::Policy::PolicyBase< dim, spacedim > > policy
::internal::DoFHandlerImplementation::NumberCache number_cache
std::vector<::internal::DoFHandlerImplementation::NumberCachemg_number_cache
std::vector< std::array< std::vector< types::global_dof_index >, dim+1 > > object_dof_indices
std::vector< std::array< std::vector< offset_type >, dim+1 > > object_dof_ptr
std::array< std::vector< types::fe_index >, dim+1 > hp_object_fe_indices
std::array< std::vector< offset_type >, dim+1 > hp_object_fe_ptr
std::vector< std::vector< types::fe_index > > hp_cell_active_fe_indices
std::vector< std::vector< types::fe_index > > hp_cell_future_fe_indices
std::vector< MGVertexDoFsmg_vertex_dofs
std::vector< std::unique_ptr<::internal::DoFHandlerImplementation::DoFLevel< dim > > > mg_levels
std::unique_ptr<::internal::DoFHandlerImplementation::DoFFaces< dim > > mg_faces
std::unique_ptr< ActiveFEIndexTransferactive_fe_index_transfer
std::vector< boost::signals2::connection > tria_listeners
std::vector< boost::signals2::connection > tria_listeners_for_transfer

Typedef Documentation

◆ level_cell_accessor

using level_cell_accessor = typename LevelSelector::CellAccessor

Definition at line 509 of file dof_handler.h.

◆ level_face_accessor

using level_face_accessor = typename LevelSelector::FaceAccessor

Definition at line 510 of file dof_handler.h.

◆ level_cell_iterator

Definition at line 512 of file dof_handler.h.

◆ level_face_iterator

Definition at line 513 of file dof_handler.h.

◆ offset_type

using offset_type = unsigned int

The type in which we store the offsets in the CRS data structures.

Definition at line 534 of file dof_handler.h.

Function Documentation

◆ DoFHandler() [1/3]

DoFHandler ( )

Standard constructor, not initializing any data. After constructing an object with this constructor, use reinit() to get a valid DoFHandler.

◆ DoFHandler() [2/3]

DoFHandler ( const Triangulation< dim, spacedim > & tria)
explicit

Constructor. Take tria as the triangulation to work on.

◆ DoFHandler() [3/3]

DoFHandler ( const DoFHandler & )
delete

Copy constructor. DoFHandler objects are large and expensive. They should not be copied, in particular not by accident, but rather deliberately constructed. As a consequence, this constructor is explicitly removed from the interface of this class.

◆ ~DoFHandler()

virtual ~DoFHandler ( )
overridevirtual

Destructor.

◆ operator=()

DoFHandler & operator= ( const DoFHandler & )
delete

Copy operator. DoFHandler objects are large and expensive. They should not be copied, in particular not by accident, but rather deliberately constructed. As a consequence, this operator is explicitly removed from the interface of this class.

◆ set_active_fe_indices() [1/2]

void set_active_fe_indices ( const std::vector< types::fe_index > & active_fe_indices)

For each locally owned cell, set the active finite element index to the corresponding value given in active_fe_indices.

The vector active_fe_indices needs to have as many entries as there are active cells. The FE indices must be in the order in which we iterate over active cells. Vector entries corresponding to active cells that are not locally owned are ignored.

Active FE indices will only be set for locally owned cells. Ghost and artificial cells will be ignored; no active FE index will be assigned to them. To exchange active FE indices on ghost cells, call distribute_dofs() afterwards.

◆ set_active_fe_indices() [2/2]

void set_active_fe_indices ( const std::vector< unsigned int > & active_fe_indices)

For each locally owned cell, set the active finite element index to the corresponding value given in active_fe_indices.

The vector active_fe_indices needs to have as many entries as there are active cells. The FE indices must be in the order in which we iterate over active cells. Vector entries corresponding to active cells that are not locally owned are ignored.

Active FE indices will only be set for locally owned cells. Ghost and artificial cells will be ignored; no active FE index will be assigned to them. To exchange active FE indices on ghost cells, call distribute_dofs() afterwards.

Deprecated
Use set_active_fe_indices() with the types::fe_index datatype.

◆ get_active_fe_indices() [1/2]

std::vector< types::fe_index > get_active_fe_indices ( ) const

For each locally relevant cell, extract the active finite element index and return them in the order in which we iterate over active cells.

As we do not know the active FE index on artificial cells, they are set to the invalid value numbers::invalid_fe_index.

For DoFHandler objects without hp-capabilities, the vector will consist of zeros, indicating that all cells use the same finite element. In hp-mode, the values may be different, though.

The returned vector has as many entries as there are active cells.

◆ get_active_fe_indices() [2/2]

void get_active_fe_indices ( std::vector< unsigned int > & active_fe_indices) const

For each locally relevant cell, extract the active finite element index and fill the vector active_fe_indices in the order in which we iterate over active cells. This vector is resized, if necessary.

As we do not know the active FE index on artificial cells, they are set to the invalid value numbers::invalid_fe_index.

For DoFHandler objects without hp-capabilities, the vector will consist of zeros, indicating that all cells use the same finite element. In hp-mode, the values may be different, though.

The returned vector has as many entries as there are active cells.

Deprecated
Use get_active_fe_indices() that returns the result vector.

◆ set_future_fe_indices()

void set_future_fe_indices ( const std::vector< types::fe_index > & future_fe_indices)

For each locally owned cell, set the future finite element index to the corresponding value given in future_fe_indices.

The vector future_fe_indices needs to have as many entries as there are active cells. The FE indices must be in the order in which we iterate over active cells. Vector entries corresponding to active cells that are not locally owned are ignored.

Future FE indices will only be set for locally owned cells. Ghost and artificial cells will be ignored; no future FE index will be assigned to them.

◆ get_future_fe_indices()

std::vector< types::fe_index > get_future_fe_indices ( ) const

For each locally owned cell, extract the future finite element index and return them in the order in which we iterate over active cells.

As we do not know the future FE index on ghost and artificial cells, they are set to the invalid value numbers::invalid_fe_index. The same applies to locally owned cells that have no future FE index assigned.

The returned vector has as many entries as there are active cells.

◆ reinit()

void reinit ( const Triangulation< dim, spacedim > & tria)

Assign a Triangulation to the DoFHandler.

Remove all associations with the previous Triangulation object and establish connections with the new one. All information about previous degrees of freedom will be removed. Activates hp-mode.

◆ distribute_dofs() [1/2]

void distribute_dofs ( const FiniteElement< dim, spacedim > & fe)

Go through the triangulation and "distribute" the degrees of freedom needed for the given finite element. "Distributing" degrees of freedom involves allocating memory to store the indices on all entities on which degrees of freedom can be located (e.g., vertices, edges, faces, etc.) and to then enumerate all degrees of freedom. In other words, while the mesh and the finite element object by themselves simply define a finite element space $V_h$, the process of distributing degrees of freedom makes sure that there is a basis for this space and that the shape functions of this basis are enumerated in an indexable, predictable way.

The exact order in which degrees of freedom on a mesh are ordered, i.e., the order in which basis functions of the finite element space are enumerated, is something that deal.II treats as an implementation detail. By and large, degrees of freedom are enumerated in the same order in which we traverse cells, but you should not rely on any specific numbering. In contrast, if you want a particular ordering, use the functions in namespace DoFRenumbering.

This function is first discussed in the introduction to the step-2 tutorial program.

Note
This function makes a copy of the finite element given as argument, and stores it as a member variable, similarly to the above function set_fe().

◆ distribute_dofs() [2/2]

void distribute_dofs ( const hp::FECollection< dim, spacedim > & fe)

Same as above but taking an hp::FECollection object.

◆ distribute_mg_dofs()

void distribute_mg_dofs ( )

Distribute level degrees of freedom on each level for geometric multigrid. The active DoFs need to be distributed using distribute_dofs() before calling this function.

◆ has_hp_capabilities()

bool has_hp_capabilities ( ) const

Returns whether this DoFHandler has hp-capabilities.

◆ has_level_dofs()

bool has_level_dofs ( ) const

This function returns whether this DoFHandler has DoFs distributed on each multigrid level or in other words if distribute_mg_dofs() has been called.

◆ has_active_dofs()

bool has_active_dofs ( ) const

This function returns whether this DoFHandler has active DoFs. This is equivalent to asking whether (i) distribute_dofs() has been called and (ii) the finite element for which degrees of freedom have been distributed actually has degrees of freedom (which is not the case for FE_Nothing, for example).

If this object is based on a parallel::distributed::Triangulation, then the current function returns true if any partition of the parallel DoFHandler object has any degrees of freedom. In other words, the function returns true even if the Triangulation does not own any active cells on the current MPI process, but at least one process owns cells and at least this one process has any degrees of freedom associated with it.

◆ initialize_local_block_info()

void initialize_local_block_info ( )

After distribute_dofs() with an FESystem element, the block structure of global and level vectors is stored in a BlockInfo object accessible with block_info(). This function initializes the local block structure on each cell in the same object.

◆ clear()

void clear ( )

Clear all data of this object.

◆ renumber_dofs() [1/2]

void renumber_dofs ( const std::vector< types::global_dof_index > & new_numbers)

Renumber degrees of freedom based on a list of new DoF indices for each of the degrees of freedom.

This function is called by the functions in DoFRenumbering function after computing a new ordering of the degree of freedom indices. However, it can of course also be called from user code.

  • new_number This array must have a size equal to the number of degrees of freedom owned by the current processor, i.e. the size must be equal to what n_locally_owned_dofs() returns. If only one processor participates in storing the current mesh, then this equals the total number of degrees of freedom, i.e. the result of n_dofs(). The contents of this array are the new global indices for each freedom listed in the IndexSet returned by locally_owned_dofs(). In the case of a sequential mesh this means that the array is a list of new indices for each of the degrees of freedom on the current mesh. In the case that we have a parallel::shared::Triangulation or parallel::distributed::Triangulation underlying this DoFHandler object, the array is a list of new indices for all the locally owned degrees of freedom, enumerated in the same order as the currently locally owned DoFs. In other words, assume that degree of freedom i is currently locally owned, then new_numbers[locally_owned_dofs().index_within_set(i)] returns the new global DoF index of i. Since the IndexSet of locally_owned_dofs() is complete in the sequential case, the latter convention for the content of the array reduces to the former in the case that only one processor participates in the mesh.
Note
While it follows from the above, it may be surprising to know that the number of locally owned degrees of freedom in a parallel computation is an invariant under renumbering, even if the indices associated with these locally owned degrees of freedom are not. At a fundamental level, this invariant exists because the decision whether a degree of freedom is locally owned or not has nothing to do with that degree of freedom's (old or new) index. Indeed, degrees of freedom are locally owned if they are on a locally owned cell and not on an interface between cells where the neighboring cell has a lower subdomain id. Since both of these conditions are independent of the index associated with the DoF, a locally owned degree of freedom will also be locally owned after renumbering. On the other hand, properties such as whether the set of indices of locally owned DoFs forms a contiguous range or not (i.e., whether the locally_owned_dofs() returns an IndexSet object for which IndexSet::is_contiguous() returns true) are of course affected by the exact renumbering performed here. For example, while the initial numbering of DoF indices done in distribute_dofs() yields a contiguous numbering, the renumbering performed by DoFRenumbering::component_wise() will, in general, not yield contiguous locally owned DoF indices.

◆ renumber_dofs() [2/2]

void renumber_dofs ( const unsigned int level,
const std::vector< types::global_dof_index > & new_numbers )

The same function as above, but renumber the degrees of freedom of a single level of a multigrid hierarchy.

◆ max_couplings_between_dofs()

unsigned int max_couplings_between_dofs ( ) const

Return the maximum number of degrees of freedom a degree of freedom in the given triangulation with the given finite element may couple with. This is the maximum number of entries per line in the system matrix; this information can therefore be used upon construction of the SparsityPattern object.

The returned number is not really the maximum number but an estimate based on the finite element and the maximum number of cells meeting at a vertex. The number holds for the constrained matrix as well.

The determination of the number of couplings can be done by simple picture drawing. An example can be found in the implementation of this function.

Note
This function is most often used to determine the maximal row length for sparsity patterns. Unfortunately, while the estimates returned by this function are rather accurate in 1d and 2d, they are often significantly too high in 3d, leading the SparsityPattern class to allocate much too much memory in some cases. Unless someone comes around to improving the present function for 3d, there is not very much one can do about these cases. The typical way to work around this problem is to use an intermediate compressed sparsity pattern that only allocates memory on demand. Refer to the step-2 and step-11 example programs on how to do this. The problem is also discussed in the documentation of the topic on Sparsity patterns.

◆ max_couplings_between_boundary_dofs()

unsigned int max_couplings_between_boundary_dofs ( ) const

Return the number of degrees of freedom located on the boundary another dof on the boundary can couple with.

The number is the same as for max_couplings_between_dofs() in one dimension less.

Note
The same applies to this function as to max_couplings_between_dofs() as regards the performance of this function. Think about one of the dynamic sparsity pattern classes instead (see Sparsity patterns).

◆ begin()

cell_iterator begin ( const unsigned int level = 0) const

Iterator to the first used cell on level level.

◆ begin_active()

active_cell_iterator begin_active ( const unsigned int level = 0) const

Iterator to the first active cell on level level. If the given level does not contain any active cells (i.e., all cells on this level are further refined), then this function returns end_active(level) so that loops of the kind

for (cell=dof_handler.begin_active(level);
cell!=dof_handler.end_active(level);
++cell)
{
...
}

have zero iterations, as may be expected if there are no active cells on this level.

◆ end() [1/2]

cell_iterator end ( ) const

Iterator past the end; this iterator serves for comparisons of iterators with past-the-end or before-the-beginning states.

◆ end() [2/2]

cell_iterator end ( const unsigned int level) const

Return an iterator which is the first iterator not on the given level. If level is the last level, then this returns end().

◆ end_active()

active_cell_iterator end_active ( const unsigned int level) const

Return an active iterator which is the first active iterator not on the given level. If level is the last level, then this returns end().

◆ begin_mg()

level_cell_iterator begin_mg ( const unsigned int level = 0) const

Iterator to the first used cell on level level. This returns a level_cell_iterator that returns level dofs when dof_indices() is called.

◆ end_mg() [1/2]

level_cell_iterator end_mg ( const unsigned int level) const

Iterator past the last cell on level level. This returns a level_cell_iterator that returns level dofs when dof_indices() is called.

◆ end_mg() [2/2]

level_cell_iterator end_mg ( ) const

Iterator past the end; this iterator serves for comparisons of iterators with past-the-end or before-the-beginning states.

◆ n_dofs() [1/2]

types::global_dof_index n_dofs ( ) const

Return the global number of degrees of freedom. If the current object handles all degrees of freedom itself (even if you may intend to solve your linear system in parallel, such as in step-17 or step-18), then this number equals the number of locally owned degrees of freedom since this object doesn't know anything about what you want to do with it and believes that it owns every degree of freedom it knows about.

On the other hand, if this object operates on a parallel::distributed::Triangulation object, then this function returns the global number of degrees of freedom, accumulated over all processors.

In either case, included in the returned number are those DoFs which are constrained by hanging nodes, see Constraints on degrees of freedom.

Mathematically speaking, the number returned by this function equals the dimension of the finite element space (without taking into account constraints) that corresponds to (i) the mesh on which it is defined, and (ii) the finite element that is used by the current object. It also, of course, equals the number of shape functions that span this space.

◆ n_dofs() [2/2]

types::global_dof_index n_dofs ( const unsigned int level) const

The (global) number of multilevel degrees of freedom on a given level.

If no level degrees of freedom have been assigned to this level, returns numbers::invalid_dof_index. Else returns the number of degrees of freedom on this level.

◆ n_boundary_dofs() [1/3]

types::global_dof_index n_boundary_dofs ( ) const

Return the number of locally owned degrees of freedom located on the boundary.

◆ n_boundary_dofs() [2/3]

template<typename number>
types::global_dof_index n_boundary_dofs ( const std::map< types::boundary_id, const Function< spacedim, number > * > & boundary_ids) const

Return the number of locally owned degrees of freedom located on those parts of the boundary which have a boundary indicator listed in the given set. The reason that a map rather than a set is used is the same as described in the documentation of that variant of DoFTools::make_boundary_sparsity_pattern() that takes a map.

There is, however, another overload of this function that takes a set argument (see below).

◆ n_boundary_dofs() [3/3]

types::global_dof_index n_boundary_dofs ( const std::set< types::boundary_id > & boundary_ids) const

Return the number of degrees of freedom located on those parts of the boundary which have a boundary indicator listed in the given set. The

◆ block_info()

const BlockInfo & block_info ( ) const

Access to an object informing of the block structure of the dof handler.

If an FESystem is used in distribute_dofs(), degrees of freedom naturally split into several blocks. For each base element as many blocks appear as its multiplicity.

At the end of distribute_dofs(), the number of degrees of freedom in each block is counted, and stored in a BlockInfo object, which can be accessed here. If you have previously called distribute_mg_dofs(), the same is done on each level of the multigrid hierarchy. Additionally, the block structure on each cell can be generated in this object by calling initialize_local_block_info().

◆ n_locally_owned_dofs()

types::global_dof_index n_locally_owned_dofs ( ) const

Return the number of degrees of freedom that belong to this process.

If this is a sequential DoFHandler, then the result equals that produced by n_dofs(). (Here, "sequential" means that either the whole program does not use MPI, or that it uses MPI but only uses a single MPI process, or that there are multiple MPI processes but the Triangulation on which this DoFHandler builds works only on one MPI process.) On the other hand, if we are operating on a parallel::distributed::Triangulation or parallel::shared::Triangulation, then it includes only the degrees of freedom that the current processor owns. Note that in this case this does not include all degrees of freedom that have been distributed on the current processor's image of the mesh: in particular, some of the degrees of freedom on the interface between the cells owned by this processor and cells owned by other processors may be theirs, and degrees of freedom on ghost cells are also not necessarily included.

◆ locally_owned_dofs()

const IndexSet & locally_owned_dofs ( ) const

Return an IndexSet describing the set of locally owned DoFs as a subset of 0..n_dofs(). The number of elements of this set equals n_locally_owned_dofs().

◆ locally_owned_mg_dofs()

const IndexSet & locally_owned_mg_dofs ( const unsigned int level) const

Return an IndexSet describing the set of locally owned DoFs used for the given multigrid level as a subset of 0..n_dofs(level).

◆ get_fe()

const FiniteElement< dim, spacedim > & get_fe ( const types::fe_index index = 0) const

Return a constant reference to the indexth finite element object that is used by this object.

◆ get_fe_collection()

const hp::FECollection< dim, spacedim > & get_fe_collection ( ) const

Return a constant reference to the set of finite element objects that are used by this object.

◆ get_triangulation()

const Triangulation< dim, spacedim > & get_triangulation ( ) const

Return a constant reference to the triangulation underlying this object.

◆ get_mpi_communicator()

MPI_Comm get_mpi_communicator ( ) const

Return MPI communicator used by the underlying triangulation.

◆ get_communicator()

MPI_Comm get_communicator ( ) const

Return MPI communicator used by the underlying triangulation.

Deprecated
Use get_mpi_communicator() instead.

◆ prepare_for_serialization_of_active_fe_indices()

void prepare_for_serialization_of_active_fe_indices ( )

Whenever serialization with a parallel::distributed::Triangulation as the underlying triangulation is considered, we also need to consider storing the active FE indices on all active cells as well.

This function registers that these indices are to be stored whenever the parallel::distributed::Triangulation::save() function is called on the underlying triangulation.

Note
Currently only implemented for triangulations of type parallel::distributed::Triangulation. An assertion will be triggered if a different type is registered.
See also
The documentation of parallel::distributed::SolutionTransfer has further information on serialization.

◆ deserialize_active_fe_indices()

void deserialize_active_fe_indices ( )

Whenever serialization with a parallel::distributed::Triangulation as the underlying triangulation is considered, we also need to consider storing the active FE indices on all active cells as well.

This function deserializes and distributes the previously stored active FE indices on all active cells.

Note
Currently only implemented for triangulations of type parallel::distributed::Triangulation. An assertion will be triggered if a different type is registered.
See also
The documentation of parallel::distributed::SolutionTransfer has further information on serialization.

◆ memory_consumption()

virtual std::size_t memory_consumption ( ) const
virtual

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

This function is made virtual, since a dof handler object might be accessed through a pointers to this base class, although the actual object might be a derived class.

◆ save()

template<class Archive>
void save ( Archive & ar,
const unsigned int version ) const

Write the data of this object to a stream for the purpose of serialization using the BOOST serialization library.

◆ load()

template<class Archive>
void load ( Archive & ar,
const unsigned int version )

Read the data of this object from a stream for the purpose of serialization using the BOOST serialization library.

◆ serialize()

template<class Archive>
void serialize ( Archive & archive,
const unsigned int version )

Write and read the data of this object from a stream for the purpose of serialization using the BOOST serialization library.

◆ clear_space()

void clear_space ( )
private

Free all memory used for non-multigrid data structures.

◆ clear_mg_space()

void clear_mg_space ( )
private

Free all memory used for multigrid data structures.

◆ setup_policy()

void setup_policy ( )
private

Set up DoFHandler policy.

◆ connect_to_triangulation_signals()

void connect_to_triangulation_signals ( )
private

Set up connections to signals of the underlying triangulation.

◆ create_active_fe_table()

void create_active_fe_table ( )
private

Create default tables for the active and future fe_indices.

Active indices are initialized with a zero indicator, meaning that fe[0] is going to be used by default. Future indices are initialized with an invalid indicator, meaning that no p-adaptation is scheduled by default.

This method is called upon construction and whenever the underlying triangulation gets created. This ensures that each cell has a valid active and future fe_index.

◆ update_active_fe_table()

void update_active_fe_table ( )
private

Update tables for active and future fe_indices.

Whenever the underlying triangulation changes (either by adaptation or deserialization), active and future FE index tables will be adjusted to the current structure of the triangulation. Missing values of active and future indices will be initialized with their defaults (see create_active_fe_table()).

This method is called post refinement and post deserialization. This ensures that each cell has a valid active and future fe_index.

◆ pre_transfer_action()

void pre_transfer_action ( )
private

A function that will be triggered through a triangulation signal just before the associated Triangulation or parallel::shared::Triangulation is modified.

The function that stores the active FE indices of all cells that will be refined or coarsened before the refinement happens, so that they can be set again after refinement.

◆ post_transfer_action()

void post_transfer_action ( )
private

A function that will be triggered through a triangulation signal just after the associated Triangulation or parallel::shared::Triangulation is modified.

The function that restores the active FE indices of all cells that were refined or coarsened.

◆ pre_distributed_transfer_action()

void pre_distributed_transfer_action ( )
private

A function that will be triggered through a triangulation signal just before the associated parallel::distributed::Triangulation is modified.

The function that stores all active FE indices on locally owned cells for distribution over all participating processors.

◆ post_distributed_transfer_action()

void post_distributed_transfer_action ( )
private

A function that will be triggered through a triangulation signal just after the associated parallel::distributed::Triangulation is modified.

The function that restores all active FE indices on locally owned cells that have been communicated.

Variable Documentation

◆ dimension

unsigned int dimension = dim
staticconstexpr

Make the dimension available in function templates.

Definition at line 519 of file dof_handler.h.

◆ space_dimension

unsigned int space_dimension = spacedim
staticconstexpr

Make the space dimension available in function templates.

Definition at line 524 of file dof_handler.h.

◆ default_fe_index

const types::fe_index default_fe_index = 0
static

The default index of the finite element to be used on a given cell.

Definition at line 529 of file dof_handler.h.

◆ block_info_object

BlockInfo block_info_object
private

An object containing information on the block structure.

Definition at line 1491 of file dof_handler.h.

◆ hp_capability_enabled

bool hp_capability_enabled
private

Boolean indicating whether or not the current DoFHandler has hp-capabilities.

Definition at line 1497 of file dof_handler.h.

◆ tria

ObserverPointer<const Triangulation<dim, spacedim>, DoFHandler<dim, spacedim> > tria
private

Address of the triangulation to work on.

Definition at line 1503 of file dof_handler.h.

◆ fe_collection

hp::FECollection<dim, spacedim> fe_collection
private

Store a hp::FECollection object. If only a single FiniteElement is used during initialization of this object, it contains the (one) FiniteElement.

Definition at line 1510 of file dof_handler.h.

◆ policy

std::unique_ptr<::internal::DoFHandlerImplementation::Policy:: PolicyBase<dim, spacedim> > policy
private

An object that describes how degrees of freedom should be distributed and renumbered.

Definition at line 1518 of file dof_handler.h.

◆ number_cache

A structure that contains all sorts of numbers that characterize the degrees of freedom this object works on.

For most members of this structure, there is an accessor function in this class that returns its value.

Definition at line 1527 of file dof_handler.h.

◆ mg_number_cache

std::vector<::internal::DoFHandlerImplementation::NumberCache> mg_number_cache
private

Data structure like number_cache, but for each multigrid level.

Definition at line 1533 of file dof_handler.h.

◆ object_dof_indices

std::vector<std::array<std::vector<types::global_dof_index>, dim + 1> > object_dof_indices
mutableprivate

Indices of degree of freedom of each d+1 geometric object (3d: vertex, line, quad, hex) for all relevant active finite elements. Identification of the appropriate position is done via object_dof_ptr (CRS scheme).

Definition at line 1541 of file dof_handler.h.

◆ object_dof_ptr

std::vector<std::array<std::vector<offset_type>, dim + 1> > object_dof_ptr
mutableprivate

Pointer to the first cached degree of freedom of a geometric object for all relevant active finite elements.

Note
In normal mode it is possible to access this data structure directly. In hp-mode, an indirection via hp_object_fe_indices/hp_object_fe_ptr is necessary.

Definition at line 1552 of file dof_handler.h.

◆ hp_object_fe_indices

std::array<std::vector<types::fe_index>, dim + 1> hp_object_fe_indices
mutableprivate

Active FE indices of each geometric object. Identification of the appropriate position of a cell in the vectors is done via hp_object_fe_ptr (CRS scheme).

Definition at line 1560 of file dof_handler.h.

◆ hp_object_fe_ptr

std::array<std::vector<offset_type>, dim + 1> hp_object_fe_ptr
mutableprivate

Pointer to the first FE index of a geometric object.

Definition at line 1565 of file dof_handler.h.

◆ hp_cell_active_fe_indices

std::vector<std::vector<types::fe_index> > hp_cell_active_fe_indices
mutableprivate

Active FE index of an active cell (identified by level and level index). This vector is only used in hp-mode.

Definition at line 1571 of file dof_handler.h.

◆ hp_cell_future_fe_indices

std::vector<std::vector<types::fe_index> > hp_cell_future_fe_indices
mutableprivate

Future FE index of an active cell (identified by level and level index). This vector is only used in hp-mode.

Definition at line 1577 of file dof_handler.h.

◆ mg_vertex_dofs

std::vector<MGVertexDoFs> mg_vertex_dofs
private

An array to store the indices for level degrees of freedom located at vertices.

Definition at line 1583 of file dof_handler.h.

◆ mg_levels

std::vector< std::unique_ptr<::internal::DoFHandlerImplementation::DoFLevel<dim> > > mg_levels
private

Space to store the DoF numbers for the different multigrid levels.

Definition at line 1590 of file dof_handler.h.

◆ mg_faces

std::unique_ptr<::internal::DoFHandlerImplementation::DoFFaces<dim> > mg_faces
private

Space to store DoF numbers of faces in the multigrid context.

Definition at line 1596 of file dof_handler.h.

◆ active_fe_index_transfer

std::unique_ptr<ActiveFEIndexTransfer> active_fe_index_transfer
private

We embed our data structure into a pointer to control that all transfer related data only exists during the actual transfer process.

Definition at line 1602 of file dof_handler.h.

◆ tria_listeners

std::vector<boost::signals2::connection> tria_listeners
private

A list of connections with which this object connects to the triangulation to get information about when the triangulation changes.

Definition at line 1608 of file dof_handler.h.

◆ tria_listeners_for_transfer

std::vector<boost::signals2::connection> tria_listeners_for_transfer
private

A list of connections with which this object connects to the triangulation. They get triggered specifically when data needs to be transferred due to refinement or repartitioning. Only active in hp-mode.

Definition at line 1615 of file dof_handler.h.