Reference documentation for deal.II version GIT relicensing-136-gb80d0be4af 2024-03-18 08:20:02+00:00
\(\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
Classes | Public Member Functions | Static Public Member Functions | Private Types | Private Member Functions | Private Attributes | Static Private Attributes | List of all members
hp::FECollection< dim, spacedim > Class Template Reference

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

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

Classes

struct  DefaultHierarchy
 

Public Member Functions

 FECollection ()
 
 FECollection (const FiniteElement< dim, spacedim > &fe)
 
template<class... FETypes>
 FECollection (const FETypes &...fes)
 
 FECollection (const std::vector< const FiniteElement< dim, spacedim > * > &fes)
 
 FECollection (const FECollection< dim, spacedim > &)=default
 
 FECollection (FECollection< dim, spacedim > &&) noexcept(std::is_nothrow_move_constructible< std::vector< std::shared_ptr< const FiniteElement< dim, spacedim > > > >::value &&std::is_nothrow_move_constructible< std::function< unsigned int(const typename hp::FECollection< dim, spacedim > &, const unsigned int)> >::value)=default
 
FECollection< dim, spacedim > & operator= (FECollection< dim, spacedim > &&)=default
 
bool operator== (const FECollection< dim, spacedim > &fe_collection) const
 
bool operator!= (const FECollection< dim, spacedim > &fe_collection) const
 
void push_back (const FiniteElement< dim, spacedim > &new_fe)
 
void push_back (const std::shared_ptr< const T > &new_entry)
 
const T & operator[] (const unsigned int index) const
 
unsigned int size () const
 
std::size_t memory_consumption () const
 
CollectionIterator< T > begin () const
 
CollectionIterator< T > end () const
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 
Querying information about the elements in the collection
unsigned int n_components () const
 
unsigned int n_blocks () const
 
unsigned int max_degree () const
 
unsigned int max_dofs_per_vertex () const
 
unsigned int max_dofs_per_line () const
 
unsigned int max_dofs_per_quad () const
 
unsigned int max_dofs_per_hex () const
 
unsigned int max_dofs_per_face () const
 
unsigned int max_dofs_per_cell () const
 
const MappingCollection< dim, spacedim > & get_reference_cell_default_linear_mapping () const
 
Functions to support hp-adaptivity
bool hp_constraints_are_implemented () const
 
std::vector< std::map< unsigned int, unsigned int > > hp_vertex_dof_identities (const std::set< unsigned int > &fes) const
 
std::vector< std::map< unsigned int, unsigned int > > hp_line_dof_identities (const std::set< unsigned int > &fes) const
 
std::vector< std::map< unsigned int, unsigned int > > hp_quad_dof_identities (const std::set< unsigned int > &fes, const unsigned int face_no=0) const
 
std::set< unsigned intfind_common_fes (const std::set< unsigned int > &fes, const unsigned int codim=0) const
 
std::set< unsigned intfind_enclosing_fes (const std::set< unsigned int > &fes, const unsigned int codim=0) const
 
unsigned int find_dominating_fe (const std::set< unsigned int > &fes, const unsigned int codim=0) const
 
unsigned int find_dominated_fe (const std::set< unsigned int > &fes, const unsigned int codim=0) const
 
unsigned int find_dominating_fe_extended (const std::set< unsigned int > &fes, const unsigned int codim=0) const
 
unsigned int find_dominated_fe_extended (const std::set< unsigned int > &fes, const unsigned int codim=0) const
 
Describing hierarchical relationships between elements
void set_hierarchy (const std::function< unsigned int(const typename hp::FECollection< dim, spacedim > &, const unsigned int)> &next, const std::function< unsigned int(const typename hp::FECollection< dim, spacedim > &, const unsigned int)> &prev)
 
void set_default_hierarchy ()
 
std::vector< unsigned intget_hierarchy_sequence (const unsigned int fe_index=0) const
 
unsigned int next_in_hierarchy (const unsigned int fe_index) const
 
unsigned int previous_in_hierarchy (const unsigned int fe_index) const
 
Components and blocks of elements
ComponentMask component_mask (const FEValuesExtractors::Scalar &scalar) const
 
ComponentMask component_mask (const FEValuesExtractors::Vector &vector) const
 
ComponentMask component_mask (const FEValuesExtractors::SymmetricTensor< 2 > &sym_tensor) const
 
ComponentMask component_mask (const BlockMask &block_mask) const
 
BlockMask block_mask (const FEValuesExtractors::Scalar &scalar) const
 
BlockMask block_mask (const FEValuesExtractors::Vector &vector) const
 
BlockMask block_mask (const FEValuesExtractors::SymmetricTensor< 2 > &sym_tensor) const
 
BlockMask block_mask (const ComponentMask &component_mask) const
 
Subscriptor functionality

Classes derived from Subscriptor provide a facility to subscribe to this object. This is mostly used by the SmartPointer class.

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
 

Static Public Member Functions

static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 
Exceptions
static ::ExceptionBaseExcNoFiniteElements ()
 

Private Types

using map_value_type = decltype(counter_map)::value_type
 
using map_iterator = decltype(counter_map)::iterator
 

Private Member Functions

void check_no_subscribers () const noexcept
 

Private Attributes

std::shared_ptr< MappingCollection< dim, spacedim > > reference_cell_default_linear_mapping
 
std::function< unsigned int(const typename hp::FECollection< dim, spacedim > &, const unsigned int)> hierarchy_next
 
std::function< unsigned int(const typename hp::FECollection< dim, spacedim > &, const unsigned int)> hierarchy_prev
 
std::vector< std::shared_ptr< const T > > entries
 
std::atomic< unsigned intcounter
 
std::map< std::string, unsigned intcounter_map
 
std::vector< std::atomic< bool > * > validity_pointers
 
const std::type_info * object_info
 

Static Private Attributes

static std::mutex mutex
 

Detailed Description

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

This class acts as a collection of finite element objects used in the DoFHandler.

It implements the concepts stated in the hp-Collections module described in the doxygen documentation.

In addition to offering access to the elements of the collection, this class provides access to the maximal number of degrees of freedom per vertex, line, etc, to allow allocation of as much memory as is necessary in the worst case when using the finite elements associated with the cells of a triangulation.

This class has not yet been implemented for the use in the codimension one case (spacedim != dim).

Definition at line 60 of file fe_collection.h.

Member Typedef Documentation

◆ map_value_type

using Subscriptor::map_value_type = decltype(counter_map)::value_type
privateinherited

The data type used in counter_map.

Definition at line 229 of file subscriptor.h.

◆ map_iterator

using Subscriptor::map_iterator = decltype(counter_map)::iterator
privateinherited

The iterator type used in counter_map.

Definition at line 234 of file subscriptor.h.

Constructor & Destructor Documentation

◆ FECollection() [1/6]

template<int dim, int spacedim>
hp::FECollection< dim, spacedim >::FECollection ( )

Default constructor. Leads to an empty collection that can later be filled using push_back(). Establishes a hierarchy of finite elements corresponding to their index in the collection.

Definition at line 31 of file fe_collection.cc.

◆ FECollection() [2/6]

template<int dim, int spacedim>
hp::FECollection< dim, spacedim >::FECollection ( const FiniteElement< dim, spacedim > &  fe)
explicit

Conversion constructor. This constructor creates a FECollection from a single finite element. More finite element objects can be added with push_back(), if desired, though it would probably be clearer to add all mappings the same way.

Definition at line 39 of file fe_collection.cc.

◆ FECollection() [3/6]

template<int dim, int spacedim>
template<class... FETypes>
hp::FECollection< dim, spacedim >::FECollection ( const FETypes &...  fes)
explicit

Constructor. This constructor creates a FECollection from one or more finite element objects passed to the constructor. For this call to be valid, all arguments need to be of types derived from class FiniteElement<dim,spacedim>.

Definition at line 880 of file fe_collection.h.

◆ FECollection() [4/6]

template<int dim, int spacedim>
hp::FECollection< dim, spacedim >::FECollection ( const std::vector< const FiniteElement< dim, spacedim > * > &  fes)

Constructor. Same as above but for any number of elements. Pointers to the elements are passed in a vector to this constructor. As above, the finite element objects pointed to by the argument are not actually used other than to create copies internally. Consequently, you can delete these pointers immediately again after calling this constructor.

Definition at line 49 of file fe_collection.cc.

◆ FECollection() [5/6]

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

Copy constructor.

◆ FECollection() [6/6]

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

Move constructor.

Note
The implementation of standard datatypes may change with different libraries, so their move members may or may not be flagged non-throwing. We need to explicitly set the noexcept specifier according to its member variables to still get the performance benefits (and to satisfy clang-tidy).

Member Function Documentation

◆ operator=()

template<int dim, int spacedim = dim>
FECollection< dim, spacedim > & hp::FECollection< dim, spacedim >::operator= ( FECollection< dim, spacedim > &&  )
default

Move assignment operator.

◆ operator==()

template<int dim, int spacedim>
bool hp::FECollection< dim, spacedim >::operator== ( const FECollection< dim, spacedim > &  fe_collection) const
inline

Equality comparison operator. All stored FiniteElement objects are compared in order.

Definition at line 919 of file fe_collection.h.

◆ operator!=()

template<int dim, int spacedim>
bool hp::FECollection< dim, spacedim >::operator!= ( const FECollection< dim, spacedim > &  fe_collection) const
inline

Non-equality comparison operator. All stored FiniteElement objects are compared in order.

Definition at line 937 of file fe_collection.h.

◆ push_back() [1/2]

template<int dim, int spacedim>
void hp::FECollection< dim, spacedim >::push_back ( const FiniteElement< dim, spacedim > &  new_fe)

Add a finite element. This function generates a copy of the given element, i.e. you can do things like push_back(FE_Q<dim>(1));. The internal copy is later destroyed by this object upon destruction of the entire collection.

When a new element is added, it needs to have the same number of vector components as all other elements already in the collection.

Definition at line 64 of file fe_collection.cc.

◆ n_components()

template<int dim, int spacedim>
unsigned int hp::FECollection< dim, spacedim >::n_components ( ) const
inline

Return the number of vector components of the finite elements in this collection. This number must be the same for all elements in the collection.

This function calls FiniteElement::n_components. See the glossary for more information.

Definition at line 900 of file fe_collection.h.

◆ n_blocks()

template<int dim, int spacedim>
unsigned int hp::FECollection< dim, spacedim >::n_blocks ( ) const

Return the number of vector blocks of the finite elements in this collection. While this class ensures that all elements stored in it have the same number of vector components, there is no such guarantees for the number of blocks each element is made up of (an element may have fewer blocks than vector components; see the glossary for more information). For example, you may have an FECollection object that stores one copy of an FESystem with dim FE_Q objects and one copy of an FE_RaviartThomas element. Both have dim vector components but while the former has dim blocks the latter has only one. Consequently, this function will throw an assertion if the number of blocks is not the same for all elements. If they are the same, this function returns the result of FiniteElement::n_blocks().

Definition at line 847 of file fe_collection.cc.

◆ max_degree()

template<int dim, int spacedim>
unsigned int hp::FECollection< dim, spacedim >::max_degree ( ) const

Return the maximum of values returned by FiniteElement::get_degree() over all elements of this collection.

Definition at line 947 of file fe_collection.h.

◆ max_dofs_per_vertex()

template<int dim, int spacedim>
unsigned int hp::FECollection< dim, spacedim >::max_dofs_per_vertex ( ) const

Return the maximal number of degrees of freedom per vertex over all elements of this collection.

Definition at line 962 of file fe_collection.h.

◆ max_dofs_per_line()

template<int dim, int spacedim>
unsigned int hp::FECollection< dim, spacedim >::max_dofs_per_line ( ) const

Return the maximal number of degrees of freedom per line over all elements of this collection.

Definition at line 977 of file fe_collection.h.

◆ max_dofs_per_quad()

template<int dim, int spacedim>
unsigned int hp::FECollection< dim, spacedim >::max_dofs_per_quad ( ) const

Return the maximal number of degrees of freedom per quad over all elements of this collection.

Definition at line 992 of file fe_collection.h.

◆ max_dofs_per_hex()

template<int dim, int spacedim>
unsigned int hp::FECollection< dim, spacedim >::max_dofs_per_hex ( ) const

Return the maximal number of degrees of freedom per hex over all elements of this collection.

Definition at line 1007 of file fe_collection.h.

◆ max_dofs_per_face()

template<int dim, int spacedim>
unsigned int hp::FECollection< dim, spacedim >::max_dofs_per_face ( ) const

Return the maximal number of degrees of freedom per face over all elements of this collection.

Definition at line 1022 of file fe_collection.h.

◆ max_dofs_per_cell()

template<int dim, int spacedim>
unsigned int hp::FECollection< dim, spacedim >::max_dofs_per_cell ( ) const

Return the maximal number of degrees of freedom per cell over all elements of this collection.

Definition at line 1037 of file fe_collection.h.

◆ get_reference_cell_default_linear_mapping()

template<int dim, int spacedim>
const MappingCollection< dim, spacedim > & hp::FECollection< dim, spacedim >::get_reference_cell_default_linear_mapping ( ) const

Return a mapping collection that consists of the default linear mappings matching the reference cells for each hp index. More details may be found in the documentation for ReferenceCell::get_default_linear_mapping().

Note
This FECollection object must remain in scope for as long as the reference cell default linear mapping is in use.

Definition at line 82 of file fe_collection.cc.

◆ hp_constraints_are_implemented()

template<int dim, int spacedim>
bool hp::FECollection< dim, spacedim >::hp_constraints_are_implemented ( ) const

Return whether all elements in this collection implement the hanging node constraints in the new way, which has to be used to make elements "hp-compatible". If this is not the case, the function returns false, which implies, that at least one element in the FECollection does not support the new face interface constraints. On the other hand, if this method does return true, this does not imply that the hp-method will work!

This behavior is related to the fact, that FiniteElement classes, which provide the new style hanging node constraints might still not provide them for all possible cases. If FE_Q and FE_RaviartThomas elements are included in the FECollection and both properly implement the get_face_interpolation_matrix method, this method will return true. But the get_face_interpolation_matrix might still fail to find an interpolation matrix between these two elements.

Definition at line 1051 of file fe_collection.h.

◆ hp_vertex_dof_identities()

template<int dim, int spacedim>
std::vector< std::map< unsigned int, unsigned int > > hp::FECollection< dim, spacedim >::hp_vertex_dof_identities ( const std::set< unsigned int > &  fes) const

This function combines the functionality of the FiniteElement::hp_vertex_dof_identities() into multi-way comparisons. Given a set of elements (whose indices are provided as argument), this function determines identities between degrees of freedom of these elements at a vertex.

The function returns a vector of such identities, where each element of the vector is a set of pairs (fe_index,dof_index) that identifies the fe_index (an element of the fes argument to this function) of an element and the dof_index indicates the how-manyth degree of freedom of that element on a vertex participates in this identity. Now, every fe_index can appear only once in these sets (for each identity, only one degree of freedom of a finite element can be involved – otherwise we would have identities between different DoFs of the same element, which would make the element not unisolvent), and as a consequence the function does not actually return a set of (fe_index,dof_index) pairs for each identity, but instead a std::map from fe_index to dof_index, which is conceptually of course equivalent to a std::set of pairs, but in practice is easier to query.

Definition at line 527 of file fe_collection.cc.

◆ hp_line_dof_identities()

template<int dim, int spacedim>
std::vector< std::map< unsigned int, unsigned int > > hp::FECollection< dim, spacedim >::hp_line_dof_identities ( const std::set< unsigned int > &  fes) const

Same as hp_vertex_dof_indices(), except that the function treats degrees of freedom on lines.

Definition at line 541 of file fe_collection.cc.

◆ hp_quad_dof_identities()

template<int dim, int spacedim>
std::vector< std::map< unsigned int, unsigned int > > hp::FECollection< dim, spacedim >::hp_quad_dof_identities ( const std::set< unsigned int > &  fes,
const unsigned int  face_no = 0 
) const

Same as hp_vertex_dof_indices(), except that the function treats degrees of freedom on quads.

Definition at line 555 of file fe_collection.cc.

◆ find_common_fes()

template<int dim, int spacedim>
std::set< unsigned int > hp::FECollection< dim, spacedim >::find_common_fes ( const std::set< unsigned int > &  fes,
const unsigned int  codim = 0 
) const

Return the indices of finite elements in this FECollection that dominate all elements associated with the provided set of indices fes.

You may find information about the domination behavior of finite elements in their respective class documentation or in the implementation of their inherited member function FiniteElement::compare_for_domination(). Consider that a finite element may or may not dominate itself (e.g. FE_Nothing elements).

For example, if a FECollection consists of {FE_Q(1),FE_Q(2),FE_Q(3),FE_Q(4)} elements and we are looking for the finite elements that dominate the middle elements of this collection (i.e., fes is {1,2}), then the answer is {FE_Q(1),FE_Q(2) and therefore this function will return their indices in the FECollection, namely {0,1}.

The codim parameter describes the codimension of the investigated subspace and specifies that it is subject to this comparison. See FiniteElement::compare_for_domination() for more information.

Definition at line 116 of file fe_collection.cc.

◆ find_enclosing_fes()

template<int dim, int spacedim>
std::set< unsigned int > hp::FECollection< dim, spacedim >::find_enclosing_fes ( const std::set< unsigned int > &  fes,
const unsigned int  codim = 0 
) const

Return the indices of finite elements in this FECollection that are dominated by all elements associated with the provided set of indices fes.

You may find information about the domination behavior of finite elements in their respective class documentation or in the implementation of their inherited member function FiniteElement::compare_for_domination(). Consider that a finite element may or may not dominate itself (e.g. FE_Nothing elements).

For example, if a FECollection consists of {FE_Q(1),FE_Q(2),FE_Q(3),FE_Q(4)} elements and we are looking for the finite elements that are dominated by the middle elements of this collection (i.e., fes is {1,2}), then the answer is {FE_Q(3),FE_Q(4) and therefore this function will return their indices in the FECollection, namely {2,3}.

The codim parameter describes the codimension of the investigated subspace and specifies that it is subject to this comparison. See FiniteElement::compare_for_domination() for more information.

Definition at line 156 of file fe_collection.cc.

◆ find_dominating_fe()

template<int dim, int spacedim>
unsigned int hp::FECollection< dim, spacedim >::find_dominating_fe ( const std::set< unsigned int > &  fes,
const unsigned int  codim = 0 
) const

Return the index of a finite element from the provided set of indices fes that dominates all other elements of this very set.

You may find information about the domination behavior of finite elements in their respective class documentation or in the implementation of their inherited member function FiniteElement::compare_for_domination(). Consider that a finite element may or may not dominate itself (e.g. FE_Nothing elements).

If this set consists of exactly one element, we consider it to be the dominating one and return its corresponding index. Further, if the function is not able to find a finite element at all, it returns numbers::invalid_fe_index.

For example, if a FECollection consists of {FE_Q(1),FE_Q(2),FE_Q(3),FE_Q(4)} elements and we are looking for the dominating finite element among the middle elements of this collection (i.e., fes is {1,2}), then the answer is FE_Q(2) and therefore this function will return its index in the FECollection, namely 1.

It is of course possible that there is more than one element that dominates all selected elements. For example, if the collection consists of {FE_Q(1),FE_Q(1),FE_Q(2),FE_Q(2)} and fes covers all indices, then one could return zero or one. In that case, the function returns either 0 or 1 since there is no tie-breaker between the two.

The codim parameter describes the codimension of the investigated subspace and specifies that it is subject to this comparison. See FiniteElement::compare_for_domination() for more information.

Definition at line 196 of file fe_collection.cc.

◆ find_dominated_fe()

template<int dim, int spacedim>
unsigned int hp::FECollection< dim, spacedim >::find_dominated_fe ( const std::set< unsigned int > &  fes,
const unsigned int  codim = 0 
) const

Return the index of a finite element from the provided set of indices fes that is dominated by all other elements of this very set.

You may find information about the domination behavior of finite elements in their respective class documentation or in the implementation of their inherited member function FiniteElement::compare_for_domination(). Consider that a finite element may or may not dominate itself (e.g. FE_Nothing elements).

If this set consists of exactly one element, we consider it to be the dominated one and return its corresponding index. Further, if the function is not able to find a finite element at all, it returns numbers::invalid_fe_index.

For example, if a FECollection consists of {FE_Q(1),FE_Q(2),FE_Q(3),FE_Q(4)} elements and we are looking for the dominated finite element among the middle elements of this collection (i.e., fes is {1,2}), then the answer is FE_Q(3) and therefore this function will return its index in the FECollection, namely 2.

It is of course possible that there is more than one element that is dominated by all selected elements. For example, if the collection consists of {FE_Q(1),FE_Q(1),FE_Q(2),FE_Q(2)} and fes covers all indices, then one could return two or three. In that case, the function returns either 2 or 3 since there is no tie-breaker between the two.

The codim parameter describes the codimension of the investigated subspace and specifies that it is subject to this comparison. See FiniteElement::compare_for_domination() for more information.

Definition at line 243 of file fe_collection.cc.

◆ find_dominating_fe_extended()

template<int dim, int spacedim>
unsigned int hp::FECollection< dim, spacedim >::find_dominating_fe_extended ( const std::set< unsigned int > &  fes,
const unsigned int  codim = 0 
) const

Return the index of a finite element from the provided set of indices fes that dominates all other elements of this very set. If we do not succeed, we extend our search on the whole collection by picking the least dominating one, which is the element that describes the largest finite element space of which all of the finite elements of the provided set fes are part of.

You may find information about the domination behavior of finite elements in their respective class documentation or in the implementation of their inherited member function FiniteElement::compare_for_domination(). Consider that a finite element may or may not dominate itself (e.g. FE_Nothing elements).

If this set consists of exactly one element, we consider it to be the dominated one and return its corresponding index. Further, if the function is not able to find a finite element at all, it returns numbers::invalid_fe_index.

The codim parameter describes the codimension of the investigated subspace and specifies that it is subject to this comparison. See FiniteElement::compare_for_domination() for more information.

Definition at line 290 of file fe_collection.cc.

◆ find_dominated_fe_extended()

template<int dim, int spacedim>
unsigned int hp::FECollection< dim, spacedim >::find_dominated_fe_extended ( const std::set< unsigned int > &  fes,
const unsigned int  codim = 0 
) const

Return the index of a finite element from the provided set of indices fes that is dominated by all other elements of this very set. If we do not succeed, we extend our search on the whole collection by picking the most dominated one, which is the element that describes the smallest finite element space which includes all finite elements of the provided set fes.

You may find information about the domination behavior of finite elements in their respective class documentation or in the implementation of their inherited member function FiniteElement::compare_for_domination(). Consider that a finite element may or may not dominate itself (e.g. FE_Nothing elements).

If this set consists of exactly one element, we consider it to be the dominating one and return its corresponding index. Further, if the function is not able to find a finite element at all, it returns numbers::invalid_fe_index.

The codim parameter describes the codimension of the investigated subspace and specifies that it is subject to this comparison. See FiniteElement::compare_for_domination() for more information.

Definition at line 310 of file fe_collection.cc.

◆ set_hierarchy()

template<int dim, int spacedim>
void hp::FECollection< dim, spacedim >::set_hierarchy ( const std::function< unsigned int(const typename hp::FECollection< dim, spacedim > &, const unsigned int)> &  next,
const std::function< unsigned int(const typename hp::FECollection< dim, spacedim > &, const unsigned int)> &  prev 
)

Set functions determining the hierarchy of finite elements, i.e. a function next that returns the index of the finite element following the given one, and a function prev returning the preceding one.

Both functions expect an hp::FECollection to be passed along with a finite element index, on whose basis the new index will be found and returned.

Note
Both passed and returned indices have to be valid within the index range of this collection, i.e. within [0, size()).

Definition at line 572 of file fe_collection.cc.

◆ set_default_hierarchy()

template<int dim, int spacedim>
void hp::FECollection< dim, spacedim >::set_default_hierarchy ( )

Set the default hierarchy corresponding to the index of each finite element in the collection.

This default hierarchy is established with functions DefaultHierarchy::next_index() and DefaultHierarchy::previous_index().

Definition at line 589 of file fe_collection.cc.

◆ get_hierarchy_sequence()

template<int dim, int spacedim>
std::vector< unsigned int > hp::FECollection< dim, spacedim >::get_hierarchy_sequence ( const unsigned int  fe_index = 0) const

Returns a sequence of FE indices that corresponds to the registered hierarchy in ascending order, i.e., FE indices are sorted from lowest to highest level.

Multiple sequences of FE indices are possible with a single custom hierarchy that can be registered with set_hierarchy(). This function will return the sequence that contains the user-provided index fe_index which could be located anywhere inside the sequence. The default hierarchy set via set_default_hierarchy(), which corresponds to FE indices in ascending order, consists of only one sequence.

This function can be used, for example, to verify that your provided hierarchy covers all elements in the desired order.

Only one sequence of FE indices exists if the size of the returned container equals the number of elements of this object, i.e., FECollection::size().

Definition at line 600 of file fe_collection.cc.

◆ next_in_hierarchy()

template<int dim, int spacedim>
unsigned int hp::FECollection< dim, spacedim >::next_in_hierarchy ( const unsigned int  fe_index) const

Function returning the index of the finite element following the given fe_index in hierarchy.

By default, the index succeeding fe_index will be returned. If fe_index already corresponds to the last index, the last index will be returned. A custom hierarchy can be supplied via the member function set_hierarchy().

Definition at line 646 of file fe_collection.cc.

◆ previous_in_hierarchy()

template<int dim, int spacedim>
unsigned int hp::FECollection< dim, spacedim >::previous_in_hierarchy ( const unsigned int  fe_index) const

Function returning the index of the finite element preceding the given fe_index in hierarchy.

By default, the index preceding fe_index will be returned. If fe_index already corresponds to the first index, the first index will be returned. A custom hierarchy can be supplied via the member function set_hierarchy().

Definition at line 661 of file fe_collection.cc.

◆ component_mask() [1/4]

template<int dim, int spacedim>
ComponentMask hp::FECollection< dim, spacedim >::component_mask ( const FEValuesExtractors::Scalar scalar) const

Return a component mask with as many elements as this object has vector components and of which exactly the one component is true that corresponds to the given argument.

Note
This function is the equivalent of FiniteElement::component_mask() with the same arguments. It verifies that it gets the same result from every one of the elements that are stored in this FECollection. If this is not the case, it throws an exception.
Parameters
scalarAn object that represents a single scalar vector component of this finite element.
Returns
A component mask that is false in all components except for the one that corresponds to the argument.

Definition at line 676 of file fe_collection.cc.

◆ component_mask() [2/4]

template<int dim, int spacedim>
ComponentMask hp::FECollection< dim, spacedim >::component_mask ( const FEValuesExtractors::Vector vector) const

Return a component mask with as many elements as this object has vector components and of which exactly the dim components are true that correspond to the given argument.

Note
This function is the equivalent of FiniteElement::component_mask() with the same arguments. It verifies that it gets the same result from every one of the elements that are stored in this FECollection. If this is not the case, it throws an exception.
Parameters
vectorAn object that represents dim vector components of this finite element.
Returns
A component mask that is false in all components except for the ones that corresponds to the argument.

Definition at line 696 of file fe_collection.cc.

◆ component_mask() [3/4]

template<int dim, int spacedim>
ComponentMask hp::FECollection< dim, spacedim >::component_mask ( const FEValuesExtractors::SymmetricTensor< 2 > &  sym_tensor) const

Return a component mask with as many elements as this object has vector components and of which exactly the dim*(dim+1)/2 components are true that correspond to the given argument.

Note
This function is the equivalent of FiniteElement::component_mask() with the same arguments. It verifies that it gets the same result from every one of the elements that are stored in this FECollection. If this is not the case, it throws an exception.
Parameters
sym_tensorAn object that represents dim*(dim+1)/2 components of this finite element that are jointly to be interpreted as forming a symmetric tensor.
Returns
A component mask that is false in all components except for the ones that corresponds to the argument.

Definition at line 716 of file fe_collection.cc.

◆ component_mask() [4/4]

template<int dim, int spacedim>
ComponentMask hp::FECollection< dim, spacedim >::component_mask ( const BlockMask block_mask) const

Given a block mask (see this glossary entry ), produce a component mask (see this glossary entry ) that represents the components that correspond to the blocks selected in the input argument. This is essentially a conversion operator from BlockMask to ComponentMask.

Note
This function is the equivalent of FiniteElement::component_mask() with the same arguments. It verifies that it gets the same result from every one of the elements that are stored in this FECollection. If this is not the case, it throws an exception.
Parameters
block_maskThe mask that selects individual blocks of the finite element
Returns
A mask that selects those components corresponding to the selected blocks of the input argument.

Definition at line 736 of file fe_collection.cc.

◆ block_mask() [1/4]

template<int dim, int spacedim>
BlockMask hp::FECollection< dim, spacedim >::block_mask ( const FEValuesExtractors::Scalar scalar) const

Return a block mask with as many elements as this object has blocks and of which exactly the one component is true that corresponds to the given argument. See the glossary for more information.

Note
This function will only succeed if the scalar referenced by the argument encompasses a complete block. In other words, if, for example, you pass an extractor for the single \(x\) velocity and this object represents an FE_RaviartThomas object, then the single scalar object you selected is part of a larger block and consequently there is no block mask that would represent it. The function will then produce an exception.
This function is the equivalent of FiniteElement::component_mask() with the same arguments. It verifies that it gets the same result from every one of the elements that are stored in this FECollection. If this is not the case, it throws an exception.
Parameters
scalarAn object that represents a single scalar vector component of this finite element.
Returns
A component mask that is false in all components except for the one that corresponds to the argument.

Definition at line 757 of file fe_collection.cc.

◆ block_mask() [2/4]

template<int dim, int spacedim>
BlockMask hp::FECollection< dim, spacedim >::block_mask ( const FEValuesExtractors::Vector vector) const

Return a component mask with as many elements as this object has vector components and of which exactly the dim components are true that correspond to the given argument. See the glossary for more information.

Note
This function is the equivalent of FiniteElement::component_mask() with the same arguments. It verifies that it gets the same result from every one of the elements that are stored in this FECollection. If this is not the case, it throws an exception.
The same caveat applies as to the version of the function above: The extractor object passed as argument must be so that it corresponds to full blocks and does not split blocks of this element.
Parameters
vectorAn object that represents dim vector components of this finite element.
Returns
A component mask that is false in all components except for the ones that corresponds to the argument.

Definition at line 779 of file fe_collection.cc.

◆ block_mask() [3/4]

template<int dim, int spacedim>
BlockMask hp::FECollection< dim, spacedim >::block_mask ( const FEValuesExtractors::SymmetricTensor< 2 > &  sym_tensor) const

Return a component mask with as many elements as this object has vector components and of which exactly the dim*(dim+1)/2 components are true that correspond to the given argument. See the glossary for more information.

Note
The same caveat applies as to the version of the function above: The extractor object passed as argument must be so that it corresponds to full blocks and does not split blocks of this element.
This function is the equivalent of FiniteElement::component_mask() with the same arguments. It verifies that it gets the same result from every one of the elements that are stored in this FECollection. If this is not the case, it throws an exception.
Parameters
sym_tensorAn object that represents dim*(dim+1)/2 components of this finite element that are jointly to be interpreted as forming a symmetric tensor.
Returns
A component mask that is false in all components except for the ones that corresponds to the argument.

Definition at line 801 of file fe_collection.cc.

◆ block_mask() [4/4]

template<int dim, int spacedim>
BlockMask hp::FECollection< dim, spacedim >::block_mask ( const ComponentMask component_mask) const

Given a component mask (see this glossary entry ), produce a block mask (see this glossary entry ) that represents the blocks that correspond to the components selected in the input argument. This is essentially a conversion operator from ComponentMask to BlockMask.

Note
This function will only succeed if the components referenced by the argument encompasses complete blocks. In other words, if, for example, you pass an component mask for the single \(x\) velocity and this object represents an FE_RaviartThomas object, then the single component you selected is part of a larger block and consequently there is no block mask that would represent it. The function will then produce an exception.
This function is the equivalent of FiniteElement::component_mask() with the same arguments. It verifies that it gets the same result from every one of the elements that are stored in this FECollection. If this is not the case, it throws an exception.
Parameters
component_maskThe mask that selects individual components of the finite element
Returns
A mask that selects those blocks corresponding to the selected blocks of the input argument.

Definition at line 824 of file fe_collection.cc.

◆ push_back() [2/2]

template<typename T >
void hp::Collection< T >::push_back ( const std::shared_ptr< const T > &  new_entry)
inherited

Add a new object.

Definition at line 255 of file collection.h.

◆ operator[]()

template<typename T >
const T & hp::Collection< T >::operator[] ( const unsigned int  index) const
inlineinherited

Return the object which was specified by the user for the active FE index which is provided as a parameter to this method.

Precondition
index must be between zero and the number of elements of the collection.

Definition at line 273 of file collection.h.

◆ size()

template<typename T >
unsigned int hp::Collection< T >::size ( ) const
inlineinherited

Return the number of objects stored in this container.

Definition at line 264 of file collection.h.

◆ memory_consumption()

template<typename T >
std::size_t hp::Collection< T >::memory_consumption ( ) const
inherited

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

Definition at line 246 of file collection.h.

◆ begin()

template<typename T >
CollectionIterator< T > hp::Collection< T >::begin ( ) const
inherited
Returns
An iterator pointing to the beginning of the underlying data (const version).

Definition at line 283 of file collection.h.

◆ end()

template<typename T >
CollectionIterator< T > hp::Collection< T >::end ( ) const
inherited
Returns
An iterator pointing to the end of the underlying data (const version).

Definition at line 292 of file collection.h.

◆ subscribe()

void Subscriptor::subscribe ( std::atomic< bool > *const  validity,
const std::string &  identifier = "" 
) const
inherited

Subscribes a user of the object by storing the pointer validity. The subscriber may be identified by text supplied as identifier.

Definition at line 135 of file subscriptor.cc.

◆ unsubscribe()

void Subscriptor::unsubscribe ( std::atomic< bool > *const  validity,
const std::string &  identifier = "" 
) const
inherited

Unsubscribes a user from the object.

Note
The identifier and the validity pointer must be the same as the one supplied to subscribe().

Definition at line 155 of file subscriptor.cc.

◆ n_subscriptions()

unsigned int Subscriptor::n_subscriptions ( ) const
inlineinherited

Return the present number of subscriptions to this object. This allows to use this class for reference counted lifetime determination where the last one to unsubscribe also deletes the object.

Definition at line 300 of file subscriptor.h.

◆ list_subscribers() [1/2]

template<typename StreamType >
void Subscriptor::list_subscribers ( StreamType &  stream) const
inlineinherited

List the subscribers to the input stream.

Definition at line 317 of file subscriptor.h.

◆ list_subscribers() [2/2]

void Subscriptor::list_subscribers ( ) const
inherited

List the subscribers to deallog.

Definition at line 203 of file subscriptor.cc.

◆ serialize()

template<class Archive >
void Subscriptor::serialize ( Archive &  ar,
const unsigned int  version 
)
inlineinherited

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

This function does not actually serialize any of the member variables of this class. The reason is that what this class stores is only who subscribes to this object, but who does so at the time of storing the contents of this object does not necessarily have anything to do with who subscribes to the object when it is restored. Consequently, we do not want to overwrite the subscribers at the time of restoring, and then there is no reason to write the subscribers out in the first place.

Definition at line 309 of file subscriptor.h.

◆ check_no_subscribers()

void Subscriptor::check_no_subscribers ( ) const
privatenoexceptinherited

Check that there are no objects subscribing to this object. If this check passes then it is safe to destroy the current object. It this check fails then this function will either abort or print an error message to deallog (by using the AssertNothrow mechanism), but will not throw an exception.

Note
Since this function is just a consistency check it does nothing in release mode.
If this function is called when there is an uncaught exception then, rather than aborting, this function prints an error message to the standard error stream and returns.

Definition at line 52 of file subscriptor.cc.

Member Data Documentation

◆ reference_cell_default_linear_mapping

template<int dim, int spacedim = dim>
std::shared_ptr<MappingCollection<dim, spacedim> > hp::FECollection< dim, spacedim >::reference_cell_default_linear_mapping
private

A linear mapping collection for all reference cell types of each index of this object.

Definition at line 855 of file fe_collection.h.

◆ hierarchy_next

template<int dim, int spacedim = dim>
std::function<unsigned int(const typename hp::FECollection<dim, spacedim> &, const unsigned int)> hp::FECollection< dim, spacedim >::hierarchy_next
private

Function returning the index of the finite element following the given one in hierarchy.

Definition at line 863 of file fe_collection.h.

◆ hierarchy_prev

template<int dim, int spacedim = dim>
std::function<unsigned int(const typename hp::FECollection<dim, spacedim> &, const unsigned int)> hp::FECollection< dim, spacedim >::hierarchy_prev
private

Function returning the index of the finite element preceding the given one in hierarchy.

Definition at line 871 of file fe_collection.h.

◆ entries

template<typename T >
std::vector<std::shared_ptr<const T> > hp::Collection< T >::entries
privateinherited

The real container, which stores pointers to the different objects.

Definition at line 236 of file collection.h.

◆ counter

std::atomic<unsigned int> Subscriptor::counter
mutableprivateinherited

Store the number of objects which subscribed to this object. Initially, this number is zero, and upon destruction it shall be zero again (i.e. all objects which subscribed should have unsubscribed again).

The creator (and owner) of an object is counted in the map below if HE manages to supply identification.

We use the mutable keyword in order to allow subscription to constant objects also.

This counter may be read from and written to concurrently in multithreaded code: hence we use the std::atomic class template.

Definition at line 218 of file subscriptor.h.

◆ counter_map

std::map<std::string, unsigned int> Subscriptor::counter_map
mutableprivateinherited

In this map, we count subscriptions for each different identification string supplied to subscribe().

Definition at line 224 of file subscriptor.h.

◆ validity_pointers

std::vector<std::atomic<bool> *> Subscriptor::validity_pointers
mutableprivateinherited

In this vector, we store pointers to the validity bool in the SmartPointer objects that subscribe to this class.

Definition at line 240 of file subscriptor.h.

◆ object_info

const std::type_info* Subscriptor::object_info
mutableprivateinherited

Pointer to the typeinfo object of this object, from which we can later deduce the class name. Since this information on the derived class is neither available in the destructor, nor in the constructor, we obtain it in between and store it here.

Definition at line 248 of file subscriptor.h.

◆ mutex

std::mutex Subscriptor::mutex
staticprivateinherited

A mutex used to ensure data consistency when accessing the mutable members of this class. This lock is used in the subscribe() and unsubscribe() functions, as well as in list_subscribers().

Definition at line 271 of file subscriptor.h.


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