Reference documentation for deal.II version 8.5.1
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]

Public Types

typedef ActiveSelector::active_cell_iterator active_cell_iterator
 
typedef ActiveSelector::cell_iterator cell_iterator
 

Public Member Functions

 DoFHandler (const Triangulation< dim, spacedim > &tria)
 
virtual ~DoFHandler ()
 
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 &&)
 
virtual ~Subscriptor ()
 
Subscriptoroperator= (const Subscriptor &)
 
Subscriptoroperator= (Subscriptor &&)
 
void subscribe (const char *identifier=0) const
 
void unsubscribe (const char *identifier=0) const
 
unsigned int n_subscriptions () 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 types::global_dof_index invalid_dof_index = numbers::invalid_dof_index
 
static const unsigned int default_fe_index = numbers::invalid_unsigned_int
 

Related Functions

(Note that these are not member 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
 
SmartPointer< const hp::FECollection< dim, spacedim >, hp::DoFHandler< dim, spacedim > > finite_elements
 
std::vector<::internal::hp::DoFLevel * > levels
 
::internal::hp::DoFIndicesOnFaces< dim > * faces
 
::internal::DoFHandler::NumberCache number_cache
 
std::vector< types::global_dof_indexvertex_dofs
 
std::vector< types::global_dof_indexvertex_dofs_offsets
 
std::vector< MGVertexDoFs > mg_vertex_dofs
 
std::vector< std::vector< bool > * > has_children
 
std::vector< boost::signals2::connection > tria_listeners
 
template<int , class , bool >
class ::DoFAccessor
 
template<class , bool >
class ::DoFCellAccessor
 
struct ::internal::DoFAccessor::Implementation
 
struct ::internal::DoFCellAccessor::Implementation
 
template<int >
class ::internal::hp::DoFIndicesOnFacesOrEdges
 
struct ::internal::hp::DoFHandler::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 hp::FECollection< dim, spacedim > & get_fe () const
 
const Triangulation< dim, spacedim > & get_tria () const 1
 
const Triangulation< dim, spacedim > & get_triangulation () const
 
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)
 
static ::ExceptionBaseExcInvalidTriangulation ()
 
static ::ExceptionBaseExcNoFESelected ()
 
static ::ExceptionBaseExcRenumberingIncomplete ()
 
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)
 
 DoFHandler (const DoFHandler &)
 
DoFHandleroperator= (const DoFHandler &)
 
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 compute_vertex_dof_identities (std::vector< types::global_dof_index > &new_dof_indices) const
 
void compute_line_dof_identities (std::vector< types::global_dof_index > &new_dof_indices) const
 
void compute_quad_dof_identities (std::vector< types::global_dof_index > &new_dof_indices) const
 
void renumber_dofs_internal (const std::vector< types::global_dof_index > &new_numbers, ::internal::int2type< 0 >)
 
void renumber_dofs_internal (const std::vector< types::global_dof_index > &new_numbers, ::internal::int2type< 1 >)
 
void renumber_dofs_internal (const std::vector< types::global_dof_index > &new_numbers, ::internal::int2type< 2 >)
 
void renumber_dofs_internal (const std::vector< types::global_dof_index > &new_numbers, ::internal::int2type< 3 >)
 

Additional Inherited Members

- Static Public Member Functions inherited from Subscriptor
static ::ExceptionBaseExcInUse (int arg1, char *arg2, std::string &arg3)
 
static ::ExceptionBaseExcNoSubscriber (char *arg1, char *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:

Author
Wolfgang Bangerth, Oliver Kayser-Herold, 2003, 2004

Definition at line 32 of file block_info.h.

Constructor & Destructor Documentation

◆ DoFHandler() [1/2]

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 1664 of file dof_handler.cc.

◆ ~DoFHandler()

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

Destructor.

Definition at line 1693 of file dof_handler.cc.

◆ DoFHandler() [2/2]

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

Copy constructor. I can see no reason why someone might want to use it, so I don't provide it. Since this class has pointer members, making it private prevents the compiler to provide it's own, incorrect one if anyone chose to copy such an object.

Member Function Documentation

◆ 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 freedoms needed for the given finite element. "Distributing" degrees of freedom involved allocating memory to store the information that describes it (e.g., whether it is located on a vertex, edge, face, etc) and to sequentially 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 purpose of this function is first discussed in the introduction to the step-2 tutorial program.

Note
A pointer of the finite element given as argument is stored. Therefore, the lifetime of the finite element object shall be longer than that of this object. If you don't want this behavior, you may want to call the clear member function which also releases the lock of this object to the finite element.

Definition at line 2574 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 2535 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 2557 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 2736 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 numbers for all the dofs.

new_numbers is an array of integers with size equal to the number of dofs on the present grid. It stores the new indices after renumbering in the order of the old indices.

This function is called by the functions in DoFRenumbering function after computing the ordering of the degrees of freedom. However, you can call this function yourself, which is necessary if a user wants to implement an ordering scheme herself, for example downwind numbering.

The new_number array must have a size equal to the number of degrees of freedom. Each entry must state the new global DoF number of the degree of freedom referenced.

Definition at line 2748 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 3099 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 3109 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 1710 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 1720 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 1736 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 1747 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 1757 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.

◆ 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 returns numbers::invalid_unsigned int 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 1964 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 section on the make_boundary_sparsity_pattern function.

The type of boundary_ids equals typename FunctionMap<spacedim,number>::type.

◆ 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

Same function, but with different data type of the argument, which is here simply a list of the boundary indicators under consideration.

Definition at line 2003 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 job, then the result equals that produced by n_dofs(). On the other hand, if we are operating on a parallel::distributed::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 job, then the vector has a single element that equals the IndexSet representing the entire range [0,n_dofs()].

◆ 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 global_dof_indices().

If this is a sequential job, then the vector has a single element equal to n_dofs().

◆ get_fe()

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.

◆ get_tria()

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

Return a constant reference to the triangulation underlying this object.

Deprecated:
Use get_triangulation() instead.

◆ 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 2067 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.

◆ operator=()

template<int dim, int spacedim = dim>
DoFHandler& hp::DoFHandler< dim, spacedim >::operator= ( const DoFHandler< dim, spacedim > &  )
private

Copy operator. I can see no reason why someone might want to use it, so I don't provide it. Since this class has pointer members, making it private prevents the compiler to provide it's own, incorrect one if anyone chose to copy such an object.

◆ clear_space()

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

Free all used memory.

Definition at line 3346 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 before distribute_dofs is called. It ensures each cell has a valid active_fe_index.

Definition at line 3147 of file dof_handler.cc.

◆ pre_refinement_action()

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

Functions that will be triggered through signals whenever the triangulation is modified.

Here they are used to administrate the the active_fe_fields during the spatial refinement.

Definition at line 3191 of file dof_handler.cc.

◆ compute_vertex_dof_identities()

template<int dim, int spacedim>
void DoFHandler< dim, spacedim >::compute_vertex_dof_identities ( std::vector< types::global_dof_index > &  new_dof_indices) const
private

Compute identities between DoFs located on vertices. Called from distribute_dofs().

Definition at line 2090 of file dof_handler.cc.

◆ compute_line_dof_identities()

template<int dim, int spacedim>
void DoFHandler< dim, spacedim >::compute_line_dof_identities ( std::vector< types::global_dof_index > &  new_dof_indices) const
private

Compute identities between DoFs located on lines. Called from distribute_dofs().

Definition at line 2224 of file dof_handler.cc.

◆ compute_quad_dof_identities()

template<int dim, int spacedim>
void DoFHandler< dim, spacedim >::compute_quad_dof_identities ( std::vector< types::global_dof_index > &  new_dof_indices) const
private

Compute identities between DoFs located on quads. Called from distribute_dofs().

Definition at line 2412 of file dof_handler.cc.

◆ renumber_dofs_internal()

template<int dim, int spacedim>
void DoFHandler< dim, spacedim >::renumber_dofs_internal ( const std::vector< types::global_dof_index > &  new_numbers,
::internal::int2type< 0 >   
)
private

Renumber the objects with the given and all lower structural dimensions, i.e. renumber vertices by giving a template argument of zero to the int2type argument, lines and vertices with one, etc.

Note that in contrast to the public renumber_dofs() function, these internal functions do not ensure that the new DoFs are contiguously numbered. The function may therefore also be used to assign different DoFs the same number, for example to unify hp DoFs corresponding to different finite elements but co-located on the same entity.

Definition at line 2799 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 889 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 898 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.

◆ 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.

◆ 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.

◆ 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.

◆ fe_is_primitive()

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

Find out whether the 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.

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 214 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 219 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 = numbers::invalid_dof_index
static

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.

Definition at line 230 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 242 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
protected

Address of the triangulation to work on.

Definition at line 708 of file dof_handler.h.

◆ finite_elements

template<int dim, int spacedim = dim>
SmartPointer<const hp::FECollection<dim,spacedim>,hp::DoFHandler<dim,spacedim> > hp::DoFHandler< dim, spacedim >::finite_elements
protected

Store a pointer to the finite element set given latest for the distribution of dofs. In order to avoid destruction of the object before the lifetime of the DoF handler, we subscribe to the finite element object. To unlock the FE before the end of the lifetime of this DoF handler, use the clear() function (this clears all data of this object as well, though).

Definition at line 718 of file dof_handler.h.

◆ levels

template<int dim, int spacedim = dim>
std::vector<::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 825 of file dof_handler.h.

◆ faces

template<int dim, int spacedim = dim>
::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 831 of file dof_handler.h.

◆ number_cache

template<int dim, int spacedim = dim>
::internal::DoFHandler::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 840 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_dofs_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 854 of file dof_handler.h.

◆ vertex_dofs_offsets

template<int dim, int spacedim = dim>
std::vector<types::global_dof_index> hp::DoFHandler< dim, spacedim >::vertex_dofs_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 868 of file dof_handler.h.

◆ has_children

template<int dim, int spacedim = dim>
std::vector<std::vector<bool> *> hp::DoFHandler< dim, spacedim >::has_children
private

Array to store the information if a cell on some level has children or not. It is used by the signal slots as a persistent buffer during the refinement, i.e. from between when pre_refinement_action is called and when post_refinement_action runs.

Definition at line 878 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 884 of file dof_handler.h.


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