Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
Classes | Public Types | Public Member Functions | Static Public Attributes | List of all members
hp::DoFHandler< dim, spacedim > Class Template Reference

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

Inheritance diagram for hp::DoFHandler< dim, spacedim >:
[legend]

Classes

struct  ActiveFEIndexTransfer
 

Public Types

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
 

Public Member Functions

 DoFHandler ()
 
 DoFHandler (const Triangulation< dim, spacedim > &tria)
 
 DoFHandler (const DoFHandler &)=delete
 
virtual ~DoFHandler () override
 
DoFHandleroperator= (const DoFHandler &)=delete
 
void initialize (const Triangulation< dim, spacedim > &tria, const hp::FECollection< dim, spacedim > &fe)
 
virtual void set_fe (const hp::FECollection< dim, spacedim > &fe)
 
virtual void distribute_dofs (const hp::FECollection< dim, spacedim > &fe)
 
void set_active_fe_indices (const std::vector< unsigned int > &active_fe_indices)
 
void get_active_fe_indices (std::vector< unsigned int > &active_fe_indices) const
 
virtual void clear ()
 
void renumber_dofs (const std::vector< types::global_dof_index > &new_numbers)
 
unsigned int max_couplings_between_dofs () const
 
unsigned int max_couplings_between_boundary_dofs () const
 
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
 
- Public Member Functions inherited from Subscriptor
 Subscriptor ()
 
 Subscriptor (const Subscriptor &)
 
 Subscriptor (Subscriptor &&) noexcept
 
virtual ~Subscriptor ()
 
Subscriptoroperator= (const Subscriptor &)
 
Subscriptoroperator= (Subscriptor &&) noexcept
 
void subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
void unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
unsigned int n_subscriptions () const
 
template<typename StreamType >
void list_subscribers (StreamType &stream) const
 
void list_subscribers () const
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 

Static Public Attributes

static const unsigned int dimension = dim
 
static const unsigned int space_dimension = spacedim
 
static const bool is_hp_dof_handler = true
 
static const types::global_dof_index invalid_dof_index
 
static const unsigned int default_fe_index = numbers::invalid_unsigned_int
 

Related Functions

(Note that these are not member functions.)

Generic Functions

Functions to support code that generically uses both DoFHandler and hp::DoFHandler.

template<int dim, int spacedim>
unsigned int max_dofs_per_cell (const hp::DoFHandler< dim, spacedim > &dh)
 
template<int dim, int spacedim>
unsigned int max_dofs_per_face (const hp::DoFHandler< dim, spacedim > &dh)
 
template<int dim, int spacedim>
unsigned int max_dofs_per_vertex (const hp::DoFHandler< dim, spacedim > &dh)
 
template<int dim, int spacedim>
unsigned int n_components (const hp::DoFHandler< dim, spacedim > &dh)
 
template<int dim, int spacedim>
bool fe_is_primitive (const hp::DoFHandler< dim, spacedim > &dh)
 

Cell iterator functions returning ranges of iterators

SmartPointer< const Triangulation< dim, spacedim >, DoFHandler< dim, spacedim > > tria
 
hp::FECollection< dim, spacedim > fe_collection
 
std::unique_ptr<::internal::DoFHandlerImplementation::Policy::PolicyBase< dim, spacedim > > policy
 
std::vector< std::unique_ptr<::internal::hp::DoFLevel > > levels
 
std::unique_ptr<::internal::hp::DoFIndicesOnFaces< dim > > faces
 
::internal::DoFHandlerImplementation::NumberCache number_cache
 
std::vector<::internal::DoFHandlerImplementation::NumberCachemg_number_cache
 
std::vector< types::global_dof_indexvertex_dofs
 
std::vector< unsigned int > vertex_dof_offsets
 
std::unique_ptr< ActiveFEIndexTransferactive_fe_index_transfer
 
std::vector< boost::signals2::connection > tria_listeners
 
template<int , class , bool >
class ::DoFAccessor
 
template<class , bool >
class ::DoFCellAccessor
 
struct ::internal::DoFAccessorImplementation::Implementation
 
struct ::internal::DoFCellAccessorImplementation::Implementation
 
template<int >
class ::internal::hp::DoFIndicesOnFacesOrEdges
 
struct ::internal::hp::DoFHandlerImplementation::Implementation
 
struct ::internal::DoFHandlerImplementation::Policy::Implementation
 
IteratorRange< cell_iteratorcell_iterators () const
 
IteratorRange< active_cell_iteratoractive_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
 
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
 
types::global_dof_index n_locally_owned_dofs () const
 
const IndexSetlocally_owned_dofs () const
 
const std::vector< IndexSet > & locally_owned_dofs_per_processor () const
 
const std::vector< types::global_dof_index > & n_locally_owned_dofs_per_processor () const
 
const IndexSetlocally_owned_mg_dofs (const unsigned int level) const
 
const std::vector< IndexSet > & locally_owned_mg_dofs_per_processor (const unsigned int level) const
 
const hp::FECollection< dim, spacedim > & get_fe () const
 
const FiniteElement< dim, spacedim > & get_fe (const unsigned int index) const
 
const hp::FECollection< dim, spacedim > & get_fe_collection () const
 
const Triangulation< dim, spacedim > & get_triangulation () const
 
virtual std::size_t memory_consumption () const
 
void prepare_for_serialization_of_active_fe_indices ()
 
void deserialize_active_fe_indices ()
 
template<class Archive >
void save (Archive &ar, const unsigned int version) const
 
template<class Archive >
void load (Archive &ar, const unsigned int version)
 
static ::ExceptionBaseExcNoFESelected ()
 
static ::ExceptionBaseExcGridsDoNotMatch ()
 
static ::ExceptionBaseExcInvalidBoundaryIndicator ()
 
static ::ExceptionBaseExcMatrixHasWrongSize (int arg1)
 
static ::ExceptionBaseExcFunctionNotUseful ()
 
static ::ExceptionBaseExcNewNumbersNotConsecutive (types::global_dof_index arg1)
 
static ::ExceptionBaseExcInvalidFEIndex (int arg1, int arg2)
 
static ::ExceptionBaseExcInvalidLevel (int arg1)
 
static ::ExceptionBaseExcFacesHaveNoLevel ()
 
static ::ExceptionBaseExcEmptyLevel (int arg1)
 
void setup_policy_and_listeners ()
 
void clear_space ()
 
template<int structdim>
types::global_dof_index get_dof_index (const unsigned int obj_level, const unsigned int obj_index, const unsigned int fe_index, const unsigned int local_index) const
 
template<int structdim>
void set_dof_index (const unsigned int obj_level, const unsigned int obj_index, const unsigned int fe_index, const unsigned int local_index, const types::global_dof_index global_index) const
 
void create_active_fe_table ()
 
void pre_refinement_action ()
 
void post_refinement_action ()
 
void pre_active_fe_index_transfer ()
 
void pre_distributed_active_fe_index_transfer ()
 
void post_active_fe_index_transfer ()
 
void post_distributed_active_fe_index_transfer ()
 
void post_distributed_serialization_of_active_fe_indices ()
 

Additional Inherited Members

- Static Public Member Functions inherited from Subscriptor
static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 

Detailed Description

template<int dim, int spacedim = dim>
class hp::DoFHandler< dim, spacedim >

Manage the distribution and numbering of the degrees of freedom for hp- FEM algorithms. This class satisfies the MeshType concept requirements.

The purpose of this class is to allow for an enumeration of degrees of freedom in the same way as the DoFHandler class, but it allows to use a different finite element on every cell. To this end, one assigns an active_fe_index to every cell that indicates which element within a collection of finite elements (represented by an object of type hp::FECollection) is the one that lives on this cell. The class then enumerates the degree of freedom associated with these finite elements on each cell of a triangulation and, if possible, identifies degrees of freedom at the interfaces of cells if they match. If neighboring cells have degrees of freedom along the common interface that do not immediate match (for example, if you have \(Q_2\) and \(Q_3\) elements meeting at a common face), then one needs to compute constraints to ensure that the resulting finite element space on the mesh remains conforming.

The whole process of working with objects of this type is explained in step-27. Many of the algorithms this class implements are described in the hp paper.

Active FE indices and their behavior under mesh refinement

The typical workflow for using this class is to create a mesh, assign an active FE index to every active cell, calls hp::DoFHandler::distribute_dofs(), and then assemble a linear system and solve a problem on this finite element space. However, one can skip assigning active FE indices upon mesh refinement in certain circumstances. In particular, the following rules apply:

Note
Finite elements need to be assigned to each cell by calling distribute_dofs() first to make this functionality available.

Active FE indices and parallel meshes

When this class is used with either a parallel::shared::Triangulation or a parallel::distributed::Triangulation, you can only set active FE indices on cells that are locally owned, using a call such as cell->set_active_fe_index(...). On the other hand, setting the active FE index on ghost or artificial cells is not allowed.

Ghost cells do acquire the information what element is active on them, however: whenever you call hp::DoFHandler::distribute_dofs(), all processors that participate in the parallel mesh exchange information in such a way that the active FE index on ghost cells equals the active FE index that was set on that processor that owned that particular ghost cell. Consequently, one can query the active_fe_index on ghost cells, just not set it by hand.

On artificial cells, no information is available about the active_fe_index used there. That's because we don't even know whether these cells exist at all, and even if they did, the current processor does not know anything specific about them. See the glossary entry on artificial cells for more information.

During refinement and coarsening, information about the active_fe_index of each cell will be automatically transferred.

However, using a parallel::distributed::Triangulation with an hp::DoFHandler requires additional attention during serialization, since no information on active FE indices will be automatically transferred. This has to be done manually using the prepare_for_serialization_of_active_fe_indices() and deserialize_active_fe_indices() functions. The former has to be called before parallel::distributed::Triangulation::save() is invoked, and the latter needs to be run after parallel::distributed::Triangulation::load(). If further data will be attached to the triangulation via the parallel::distributed::CellDataTransfer, parallel::distributed::SolutionTransfer, or Particles::ParticleHandler classes, all corresponding preparation and deserialization function calls need to happen in the same order. Consult the documentation of parallel::distributed::SolutionTransfer for more information.

Author
Wolfgang Bangerth, 2003, 2004, 2017, 2018
Oliver Kayser-Herold, 2003, 2004
Marc Fehling, 2018

Definition at line 35 of file block_info.h.

Constructor & Destructor Documentation

◆ DoFHandler() [1/3]

template<int dim, int spacedim>
DoFHandler< dim, spacedim >::DoFHandler ( )

Default Constructor.

Definition at line 1246 of file dof_handler.cc.

◆ DoFHandler() [2/3]

template<int dim, int spacedim>
DoFHandler< dim, spacedim >::DoFHandler ( const Triangulation< dim, spacedim > &  tria)

Constructor. Take tria as the triangulation to work on.

Definition at line 1254 of file dof_handler.cc.

◆ DoFHandler() [3/3]

template<int dim, int spacedim = dim>
hp::DoFHandler< dim, spacedim >::DoFHandler ( const DoFHandler< dim, spacedim > &  )
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()

template<int dim, int spacedim>
DoFHandler< dim, spacedim >::~DoFHandler ( )
overridevirtual

Destructor.

Definition at line 1267 of file dof_handler.cc.

Member Function Documentation

◆ operator=()

template<int dim, int spacedim = dim>
DoFHandler& hp::DoFHandler< dim, spacedim >::operator= ( const DoFHandler< dim, spacedim > &  )
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.

◆ initialize()

template<int dim, int spacedim>
void DoFHandler< dim, spacedim >::initialize ( const Triangulation< dim, spacedim > &  tria,
const hp::FECollection< dim, spacedim > &  fe 
)

Assign a Triangulation and a FECollection to the DoFHandler and compute the distribution of degrees of freedom over the mesh.

Definition at line 1557 of file dof_handler.cc.

◆ set_fe()

template<int dim, int spacedim>
void DoFHandler< dim, spacedim >::set_fe ( const hp::FECollection< dim, spacedim > &  fe)
virtual

Assign a hp::FECollection fe to this object.

In case a parallel::Triangulation is assigned to this object, the active_fe_indices will be exchanged between processors so that each one knows the indices on its own cells and all ghost cells.

Note
In accordance with DoFHandler::set_fe(), this function also makes a copy of the object given as argument.
Warning
This function only sets a hp::FECollection. Degrees of freedom have either not been distributed yet, or are distributed using a previously set collection. In both cases, accessing degrees of freedom will lead to invalid results. To restore consistency, call distribute_dofs().

Definition at line 1581 of file dof_handler.cc.

◆ distribute_dofs()

template<int dim, int spacedim>
void DoFHandler< dim, spacedim >::distribute_dofs ( const hp::FECollection< dim, spacedim > &  fe)
virtual

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.

Note
In accordance with DoFHandler::distribute_dofs(), this function also makes a copy of the object given as argument.

Definition at line 1617 of file dof_handler.cc.

◆ set_active_fe_indices()

template<int dim, int spacedim>
void DoFHandler< dim, spacedim >::set_active_fe_indices ( const std::vector< unsigned int > &  active_fe_indices)

Go through the triangulation and set the active FE indices of all active cells to the values given in active_fe_indices.

Definition at line 1519 of file dof_handler.cc.

◆ get_active_fe_indices()

template<int dim, int spacedim>
void DoFHandler< dim, spacedim >::get_active_fe_indices ( std::vector< unsigned int > &  active_fe_indices) const

Go through the triangulation and store the active FE indices of all active cells to the vector active_fe_indices. This vector is resized, if necessary.

Definition at line 1540 of file dof_handler.cc.

◆ clear()

template<int dim, int spacedim>
void DoFHandler< dim, spacedim >::clear ( )
virtual

Clear all data of this object and especially delete the lock this object has to the finite element used the last time when distribute_dofs was called.

Definition at line 1793 of file dof_handler.cc.

◆ renumber_dofs()

template<int dim, int spacedim>
void DoFHandler< dim, spacedim >::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.

Definition at line 1803 of file dof_handler.cc.

◆ max_couplings_between_dofs()

template<int dim, int spacedim>
unsigned int DoFHandler< dim, spacedim >::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 also.

As for DoFHandler::max_couplings_between_dofs(), the result of this function is often not very accurate for 3d and/or high polynomial degrees. The consequences are discussed in the documentation of the module on Sparsity patterns.

Definition at line 1870 of file dof_handler.cc.

◆ max_couplings_between_boundary_dofs()

template<int dim, int spacedim>
unsigned int DoFHandler< dim, spacedim >::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_coupling_between_dofs in one dimension less.

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

Definition at line 1881 of file dof_handler.cc.

◆ begin()

template<int dim, int spacedim>
DoFHandler< dim, spacedim >::cell_iterator DoFHandler< dim, spacedim >::begin ( const unsigned int  level = 0) const

Iterator to the first used cell on level level.

Definition at line 1289 of file dof_handler.cc.

◆ begin_active()

template<int dim, int spacedim>
DoFHandler< dim, spacedim >::active_cell_iterator DoFHandler< dim, spacedim >::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.

Definition at line 1298 of file dof_handler.cc.

◆ end() [1/2]

template<int dim, int spacedim>
DoFHandler< dim, spacedim >::cell_iterator DoFHandler< dim, spacedim >::end ( ) const

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

Definition at line 1314 of file dof_handler.cc.

◆ end() [2/2]

template<int dim, int spacedim>
DoFHandler< dim, spacedim >::cell_iterator DoFHandler< dim, spacedim >::end ( const unsigned int  level) const

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

Definition at line 1323 of file dof_handler.cc.

◆ end_active()

template<int dim, int spacedim>
DoFHandler< dim, spacedim >::active_cell_iterator DoFHandler< dim, spacedim >::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().

Definition at line 1334 of file dof_handler.cc.

◆ n_dofs() [1/2]

template<int dim, int spacedim = dim>
types::global_dof_index hp::DoFHandler< dim, spacedim >::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]

template<int dim, int spacedim = dim>
types::global_dof_index hp::DoFHandler< dim, spacedim >::n_dofs ( const unsigned int  level) const

The number of multilevel dofs on given level. Since hp::DoFHandler does not support multilevel methods yet, this function throws an exception ExcNotImplemented() independent of its argument.

◆ n_boundary_dofs() [1/3]

template<int dim, int spacedim>
types::global_dof_index DoFHandler< dim, spacedim >::n_boundary_dofs ( ) const

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

Definition at line 1392 of file dof_handler.cc.

◆ n_boundary_dofs() [2/3]

template<int dim, int spacedim = dim>
template<typename number >
types::global_dof_index hp::DoFHandler< dim, spacedim >::n_boundary_dofs ( const std::map< types::boundary_id, const Function< spacedim, number > *> &  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 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]

template<int dim, int spacedim>
types::global_dof_index DoFHandler< dim, spacedim >::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

Definition at line 1427 of file dof_handler.cc.

◆ n_locally_owned_dofs()

template<int dim, int spacedim = dim>
types::global_dof_index hp::DoFHandler< dim, spacedim >::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()

template<int dim, int spacedim = dim>
const IndexSet& hp::DoFHandler< dim, spacedim >::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_dofs_per_processor()

template<int dim, int spacedim = dim>
const std::vector<IndexSet>& hp::DoFHandler< dim, spacedim >::locally_owned_dofs_per_processor ( ) const

Return a vector that stores the locally owned DoFs of each processor. If you are only interested in the number of elements each processor owns then n_dofs_per_processor() is a better choice.

If this is a sequential DoFHandler, then the vector has a single element that equals the IndexSet representing the entire range [0,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.)

◆ n_locally_owned_dofs_per_processor()

template<int dim, int spacedim = dim>
const std::vector<types::global_dof_index>& hp::DoFHandler< dim, spacedim >::n_locally_owned_dofs_per_processor ( ) const

Return a vector that stores the number of degrees of freedom each processor that participates in this triangulation owns locally. The sum of all these numbers equals the number of degrees of freedom that exist globally, i.e. what n_dofs() returns.

Each element of the vector returned by this function equals the number of elements of the corresponding sets returned by locally_owned_dofs_per_processor().

If this is a sequential DoFHandler, then the vector has a single element equal to 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.)

◆ locally_owned_mg_dofs()

template<int dim, int spacedim = dim>
const IndexSet& hp::DoFHandler< dim, spacedim >::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. Since hp::DoFHandler does not support multilevel methods yet, this function throws an exception ExcNotImplemented() independent of its argument.

◆ locally_owned_mg_dofs_per_processor()

template<int dim, int spacedim = dim>
const std::vector<IndexSet>& hp::DoFHandler< dim, spacedim >::locally_owned_mg_dofs_per_processor ( const unsigned int  level) const

Return a vector that stores the locally owned level DoFs of each processor on the given level level. Since hp::DoFHandler does not support multilevel methods yet, this function throws an exception ExcNotImplemented() independent of its argument.

◆ get_fe() [1/2]

template<int dim, int spacedim = dim>
const hp::FECollection<dim, spacedim>& hp::DoFHandler< dim, spacedim >::get_fe ( ) const

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

Deprecated:
Use get_fe_collection() instead.

◆ get_fe() [2/2]

template<int dim, int spacedim = dim>
const FiniteElement<dim, spacedim>& hp::DoFHandler< dim, spacedim >::get_fe ( const unsigned int  index) const

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

◆ get_fe_collection()

template<int dim, int spacedim = dim>
const hp::FECollection<dim, spacedim>& hp::DoFHandler< dim, spacedim >::get_fe_collection ( ) const

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

◆ get_triangulation()

template<int dim, int spacedim = dim>
const Triangulation<dim, spacedim>& hp::DoFHandler< dim, spacedim >::get_triangulation ( ) const

Return a constant reference to the triangulation underlying this object.

◆ memory_consumption()

template<int dim, int spacedim>
std::size_t DoFHandler< dim, spacedim >::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.

Definition at line 1497 of file dof_handler.cc.

◆ prepare_for_serialization_of_active_fe_indices()

template<int dim, int spacedim>
void DoFHandler< dim, spacedim >::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.

Definition at line 2160 of file dof_handler.cc.

◆ deserialize_active_fe_indices()

template<int dim, int spacedim>
void DoFHandler< dim, spacedim >::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.

Definition at line 2236 of file dof_handler.cc.

◆ save()

template<int dim, int spacedim = dim>
template<class Archive >
void hp::DoFHandler< dim, spacedim >::save ( Archive &  ar,
const unsigned int  version 
) const

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

◆ load()

template<int dim, int spacedim = dim>
template<class Archive >
void hp::DoFHandler< dim, spacedim >::load ( Archive &  ar,
const unsigned int  version 
)

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

◆ setup_policy_and_listeners()

template<int dim, int spacedim>
void DoFHandler< dim, spacedim >::setup_policy_and_listeners ( )
private

Setup policy and listeners based on the underlying Triangulation.

Definition at line 1692 of file dof_handler.cc.

◆ clear_space()

template<int dim, int spacedim>
void DoFHandler< dim, spacedim >::clear_space ( )
private

Free all used memory.

Definition at line 2322 of file dof_handler.cc.

◆ create_active_fe_table()

template<int dim, int spacedim>
void DoFHandler< dim, spacedim >::create_active_fe_table ( )
private

Create default tables for the active_fe_indices in the internal::hp::DoFLevel. They are initialized with a zero indicator, meaning that fe[0] is going to be used by default. This method is called before refinement and while setting the finite elements via set_fe(). It ensures each cell has a valid active_fe_index.

Definition at line 1915 of file dof_handler.cc.

◆ pre_refinement_action()

template<int dim, int spacedim>
void DoFHandler< dim, spacedim >::pre_refinement_action ( )
private

A function that will be triggered through a triangulation signal just before the triangulation is modified.

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

Definition at line 1960 of file dof_handler.cc.

◆ post_refinement_action()

template<int dim, int spacedim>
void DoFHandler< dim, spacedim >::post_refinement_action ( )
private

A function that will be triggered through a triangulation signal just after the triangulation is modified.

The function that restores the active_fe_flags of all cells that were refined.

Definition at line 1969 of file dof_handler.cc.

◆ pre_active_fe_index_transfer()

template<int dim, int spacedim>
void DoFHandler< dim, spacedim >::pre_active_fe_index_transfer ( )
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.

Definition at line 2004 of file dof_handler.cc.

◆ pre_distributed_active_fe_index_transfer()

template<int dim, int spacedim>
void DoFHandler< dim, spacedim >::pre_distributed_active_fe_index_transfer ( )
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.

Definition at line 2024 of file dof_handler.cc.

◆ post_active_fe_index_transfer()

template<int dim, int spacedim>
void DoFHandler< dim, spacedim >::post_active_fe_index_transfer ( )
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.

Definition at line 2104 of file dof_handler.cc.

◆ post_distributed_active_fe_index_transfer()

template<int dim, int spacedim>
void DoFHandler< dim, spacedim >::post_distributed_active_fe_index_transfer ( )
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.

Definition at line 2130 of file dof_handler.cc.

◆ post_distributed_serialization_of_active_fe_indices()

template<int dim, int spacedim>
void DoFHandler< dim, spacedim >::post_distributed_serialization_of_active_fe_indices ( )
private

A function that will be triggered through a triangulation signal just after the associated parallel::distributed::Triangulation has been saved.

The function frees all memory related to the transfer of active_fe_indices.

Definition at line 2216 of file dof_handler.cc.

Friends And Related Function Documentation

◆ ::DoFAccessor

template<int dim, int spacedim = dim>
template<int , class , bool >
friend class ::DoFAccessor
friend

Make accessor objects friends.

Definition at line 1289 of file dof_handler.h.

◆ ::internal::hp::DoFIndicesOnFacesOrEdges

template<int dim, int spacedim = dim>
template<int >
friend class ::internal::hp::DoFIndicesOnFacesOrEdges
friend

Likewise for DoFLevel objects since they need to access the vertex dofs in the functions that set and retrieve vertex dof indices.

Definition at line 1301 of file dof_handler.h.

◆ max_dofs_per_cell()

template<int dim, int spacedim>
unsigned int max_dofs_per_cell ( const hp::DoFHandler< dim, spacedim > &  dh)
related

Maximal number of degrees of freedom on a cell.

Deprecated:
Use dh.get_fe_collection().max_dofs_per_cell().

◆ max_dofs_per_face()

template<int dim, int spacedim>
unsigned int max_dofs_per_face ( const hp::DoFHandler< dim, spacedim > &  dh)
related

Maximal number of degrees of freedom on a face.

This function exists for both non-hp and hp DoFHandlers, to allow for a uniform interface to query this property.

Deprecated:
Use dh.get_fe_collection().max_dofs_per_face().

◆ max_dofs_per_vertex()

template<int dim, int spacedim>
unsigned int max_dofs_per_vertex ( const hp::DoFHandler< dim, spacedim > &  dh)
related

Maximal number of degrees of freedom on a vertex.

This function exists for both non-hp and hp DoFHandlers, to allow for a uniform interface to query this property.

Deprecated:
Use dh.get_fe_collection().max_dofs_per_vertex().

◆ n_components()

template<int dim, int spacedim>
unsigned int n_components ( const hp::DoFHandler< dim, spacedim > &  dh)
related

Number of vector components in the finite element object used by this DoFHandler.

This function exists for both non-hp and hp DoFHandlers, to allow for a uniform interface to query this property.

Deprecated:
Use dh.get_fe_collection().n_components().

◆ fe_is_primitive()

template<int dim, int spacedim>
bool fe_is_primitive ( const hp::DoFHandler< dim, spacedim > &  dh)
related

Find out whether the first FiniteElement used by this DoFHandler is primitive or not.

This function exists for both non-hp and hp DoFHandlers, to allow for a uniform interface to query this property.

Deprecated:
Use dh.get_fe(0).is_primitive().

Member Data Documentation

◆ dimension

template<int dim, int spacedim = dim>
const unsigned int DoFHandler< dim, spacedim >::dimension = dim
static

Make the dimension available in function templates.

Definition at line 342 of file dof_handler.h.

◆ space_dimension

template<int dim, int spacedim = dim>
const unsigned int hp::DoFHandler< dim, spacedim >::space_dimension = spacedim
static

Make the space dimension available in function templates.

Definition at line 347 of file dof_handler.h.

◆ is_hp_dof_handler

template<int dim, int spacedim = dim>
const bool hp::DoFHandler< dim, spacedim >::is_hp_dof_handler = true
static

Make the type of this DoFHandler available in function templates.

Definition at line 352 of file dof_handler.h.

◆ invalid_dof_index

template<int dim, int spacedim = dim>
const types::global_dof_index DoFHandler< dim, spacedim >::invalid_dof_index
static
Initial value:

When the arrays holding the DoF indices are set up, but before they are filled with actual values, they are set to an invalid value, in order to monitor possible problems. This invalid value is the constant defined here.

Please note that you should not rely on it having a certain value, but rather take its symbolic name.

Deprecated:
Use numbers::invalid_dof_index instead.

Definition at line 366 of file dof_handler.h.

◆ default_fe_index

template<int dim, int spacedim = dim>
const unsigned int DoFHandler< dim, spacedim >::default_fe_index = numbers::invalid_unsigned_int
static

The default index of the finite element to be used on a given cell. For the usual, non-hp DoFHandler class that only supports the same finite element to be used on all cells, the index of the finite element needs to be the same on all cells anyway, and by convention we pick zero for this value. The situation here is different, since the hp classes support the case where different finite element indices may be used on different cells. The default index consequently corresponds to an invalid value.

Definition at line 379 of file dof_handler.h.

◆ tria

template<int dim, int spacedim = dim>
SmartPointer<const Triangulation<dim, spacedim>, DoFHandler<dim, spacedim> > hp::DoFHandler< dim, spacedim >::tria
private

Address of the triangulation to work on.

Definition at line 1040 of file dof_handler.h.

◆ fe_collection

template<int dim, int spacedim = dim>
hp::FECollection<dim, spacedim> hp::DoFHandler< dim, spacedim >::fe_collection
private

Store a copy of the finite element set given latest to distribute_dofs().

Definition at line 1045 of file dof_handler.h.

◆ policy

template<int dim, int spacedim = dim>
std::unique_ptr<::internal::DoFHandlerImplementation::Policy:: PolicyBase<dim, spacedim> > hp::DoFHandler< dim, spacedim >::policy
private

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

Definition at line 1053 of file dof_handler.h.

◆ levels

template<int dim, int spacedim = dim>
std::vector<std::unique_ptr<::internal::hp::DoFLevel> > hp::DoFHandler< dim, spacedim >::levels
private

Space to store the DoF numbers for the different levels. Analogous to the levels[] tree of the Triangulation objects.

Definition at line 1174 of file dof_handler.h.

◆ faces

template<int dim, int spacedim = dim>
std::unique_ptr<::internal::hp::DoFIndicesOnFaces<dim> > hp::DoFHandler< dim, spacedim >::faces
private

Space to store the DoF numbers for the faces. Analogous to the faces pointer of the Triangulation objects.

Definition at line 1180 of file dof_handler.h.

◆ number_cache

template<int dim, int spacedim = dim>
::internal::DoFHandlerImplementation::NumberCache hp::DoFHandler< dim, spacedim >::number_cache
private

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 1189 of file dof_handler.h.

◆ mg_number_cache

template<int dim, int spacedim = dim>
std::vector<::internal::DoFHandlerImplementation::NumberCache> hp::DoFHandler< dim, spacedim >::mg_number_cache
private

A structure that contains all sorts of numbers that characterize the degrees of freedom on multigrid levels. Since multigrid is not currently supported, this table is not filled with valid entries.

Definition at line 1197 of file dof_handler.h.

◆ vertex_dofs

template<int dim, int spacedim = dim>
std::vector<types::global_dof_index> hp::DoFHandler< dim, spacedim >::vertex_dofs
private

Array to store the indices for degrees of freedom located at vertices.

The format used here, in the form of a linked list, is the same as used for the arrays used in the internal::hp::DoFLevel hierarchy. Starting indices into this array are provided by the vertex_dof_offsets field.

Access to this field is generally through the DoFAccessor::get_vertex_dof_index() and DoFAccessor::set_vertex_dof_index() functions, encapsulating the actual data format used to the present class.

Definition at line 1211 of file dof_handler.h.

◆ vertex_dof_offsets

template<int dim, int spacedim = dim>
std::vector<unsigned int> hp::DoFHandler< dim, spacedim >::vertex_dof_offsets
private

For each vertex in the triangulation, store the offset within the vertex_dofs array where the dofs for this vertex start.

As for that array, the format is the same as described in the documentation of hp::DoFLevel.

Access to this field is generally through the Accessor::get_vertex_dof_index() and Accessor::set_vertex_dof_index() functions, encapsulating the actual data format used to the present class.

Definition at line 1225 of file dof_handler.h.

◆ active_fe_index_transfer

template<int dim, int spacedim = dim>
std::unique_ptr<ActiveFEIndexTransfer> hp::DoFHandler< dim, spacedim >::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 1277 of file dof_handler.h.

◆ tria_listeners

template<int dim, int spacedim = dim>
std::vector<boost::signals2::connection> hp::DoFHandler< dim, spacedim >::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 1283 of file dof_handler.h.


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