Reference documentation for deal.II version 9.2.0
\(\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\}}\)
Classes | Public Member Functions | List of all members

#include <deal.II/fe/fe.h>

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

Classes

class  InternalData
 

Public Member Functions

 FESystem (const FiniteElement< dim, spacedim > &fe, const unsigned int n_elements)
 
 FESystem (const FiniteElement< dim, spacedim > &fe1, const unsigned int n1, const FiniteElement< dim, spacedim > &fe2, const unsigned int n2)
 
 FESystem (const FiniteElement< dim, spacedim > &fe1, const unsigned int n1, const FiniteElement< dim, spacedim > &fe2, const unsigned int n2, const FiniteElement< dim, spacedim > &fe3, const unsigned int n3)
 
 FESystem (const FiniteElement< dim, spacedim > &fe1, const unsigned int n1, const FiniteElement< dim, spacedim > &fe2, const unsigned int n2, const FiniteElement< dim, spacedim > &fe3, const unsigned int n3, const FiniteElement< dim, spacedim > &fe4, const unsigned int n4)
 
 FESystem (const FiniteElement< dim, spacedim > &fe1, const unsigned int n1, const FiniteElement< dim, spacedim > &fe2, const unsigned int n2, const FiniteElement< dim, spacedim > &fe3, const unsigned int n3, const FiniteElement< dim, spacedim > &fe4, const unsigned int n4, const FiniteElement< dim, spacedim > &fe5, const unsigned int n5)
 
 FESystem (const std::vector< const FiniteElement< dim, spacedim > * > &fes, const std::vector< unsigned int > &multiplicities)
 
template<class... FEPairs, typename = typename enable_if_all< (std::is_same<typename std::decay<FEPairs>::type, std::pair<std::unique_ptr<FiniteElement<dim, spacedim>>, unsigned int>>::value || std::is_base_of<FiniteElement<dim, spacedim>, typename std::decay<FEPairs>::type>::value)...>::type>
 FESystem (FEPairs &&... fe_pairs)
 
 FESystem (const std::initializer_list< std::pair< std::unique_ptr< FiniteElement< dim, spacedim >>, unsigned int >> &fe_systems)
 
 FESystem (const FESystem< dim, spacedim > &)=delete
 
 FESystem (FESystem< dim, spacedim > &&other_fe_system) noexcept
 
virtual ~FESystem () override=default
 
virtual std::string get_name () const override
 
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone () const override
 
virtual UpdateFlags requires_update_flags (const UpdateFlags update_flags) const override
 
virtual const FiniteElement< dim, spacedim > & get_sub_fe (const unsigned int first_component, const unsigned int n_selected_components) const override
 
virtual double shape_value (const unsigned int i, const Point< dim > &p) const override
 
virtual double shape_value_component (const unsigned int i, const Point< dim > &p, const unsigned int component) const override
 
virtual Tensor< 1, dim > shape_grad (const unsigned int i, const Point< dim > &p) const override
 
virtual Tensor< 1, dim > shape_grad_component (const unsigned int i, const Point< dim > &p, const unsigned int component) const override
 
virtual Tensor< 2, dim > shape_grad_grad (const unsigned int i, const Point< dim > &p) const override
 
virtual Tensor< 2, dim > shape_grad_grad_component (const unsigned int i, const Point< dim > &p, const unsigned int component) const override
 
virtual Tensor< 3, dim > shape_3rd_derivative (const unsigned int i, const Point< dim > &p) const override
 
virtual Tensor< 3, dim > shape_3rd_derivative_component (const unsigned int i, const Point< dim > &p, const unsigned int component) const override
 
virtual Tensor< 4, dim > shape_4th_derivative (const unsigned int i, const Point< dim > &p) const override
 
virtual Tensor< 4, dim > shape_4th_derivative_component (const unsigned int i, const Point< dim > &p, const unsigned int component) const override
 
virtual void get_interpolation_matrix (const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const override
 
virtual const FiniteElement< dim, spacedim > & base_element (const unsigned int index) const override
 
virtual bool has_support_on_face (const unsigned int shape_index, const unsigned int face_index) const override
 
virtual const FullMatrix< double > & get_restriction_matrix (const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const override
 
virtual const FullMatrix< double > & get_prolongation_matrix (const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const override
 
virtual unsigned int face_to_cell_index (const unsigned int face_dof_index, const unsigned int face, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false) const override
 
virtual Point< dim > unit_support_point (const unsigned int index) const override
 
virtual Point< dim - 1 > unit_face_support_point (const unsigned int index) const override
 
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes () const override
 

Functions to support hp

static const unsigned int invalid_face_number = numbers::invalid_unsigned_int
 
std::vector< std::pair< std::unique_ptr< const FiniteElement< dim, spacedim > >, unsigned int > > base_elements
 
std::vector< std::vector< std::size_t > > generalized_support_points_index_table
 
std::mutex mutex
 
class FE_Enriched< dim, spacedim >
 
virtual bool hp_constraints_are_implemented () const override
 
virtual void get_face_interpolation_matrix (const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const override
 
virtual void get_subface_interpolation_matrix (const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix) const override
 
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities (const FiniteElement< dim, spacedim > &fe_other) const override
 
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities (const FiniteElement< dim, spacedim > &fe_other) const override
 
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_quad_dof_identities (const FiniteElement< dim, spacedim > &fe_other) const override
 
virtual FiniteElementDomination::Domination compare_for_domination (const FiniteElement< dim, spacedim > &fe_other, const unsigned int codim=0) const override final
 
virtual void convert_generalized_support_point_values_to_dof_values (const std::vector< Vector< double >> &support_point_values, std::vector< double > &dof_values) const override
 
virtual std::size_t memory_consumption () const override
 
virtual std::unique_ptr< typename FiniteElement< dim, spacedim >::InternalDataBase > get_data (const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
 
virtual std::unique_ptr< typename FiniteElement< dim, spacedim >::InternalDataBase > get_face_data (const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim - 1 > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
 
virtual std::unique_ptr< typename FiniteElement< dim, spacedim >::InternalDataBase > get_subface_data (const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim - 1 > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
 
virtual void fill_fe_values (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const Quadrature< dim > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const typename FiniteElement< dim, spacedim >::InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
 
virtual void fill_fe_face_values (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const Quadrature< dim - 1 > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const typename FiniteElement< dim, spacedim >::InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
 
virtual void fill_fe_subface_values (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int sub_no, const Quadrature< dim - 1 > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const typename FiniteElement< dim, spacedim >::InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
 
template<int dim_1>
void compute_fill (const Mapping< dim, spacedim > &mapping, const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int sub_no, const Quadrature< dim_1 > &quadrature, const CellSimilarity::Similarity cell_similarity, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const typename FiniteElement< dim, spacedim >::InternalDataBase &fe_data, const internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const
 
void initialize (const std::vector< const FiniteElement< dim, spacedim > * > &fes, const std::vector< unsigned int > &multiplicities)
 
void build_interface_constraints ()
 
template<int structdim>
std::vector< std::pair< unsigned int, unsigned int > > hp_object_dof_identities (const FiniteElement< dim, spacedim > &fe_other) const
 

Detailed Description

template<int dim, int spacedim = dim>
class FESystem< dim, spacedim >

This class provides an interface to group several elements together into one. To the outside world, the resulting object looks just like a usual finite element object, which is composed of several other finite elements that are possibly of different type. The result is then a vector-valued finite element. An example is given in the documentation of namespace FETools::Compositing, when using the "tensor product" strategy.

Vector valued elements are discussed in a number of tutorial programs, for example step-8, step-20, step-21, and in particular in the Handling vector valued problems module.

Note
The material presented here is also discussed in video lecture 19, video lecture 20. (All video lectures are also available here.)

FESystem, components and blocks

An FESystem, except in the most trivial case, produces a vector-valued finite element with several components. The number of components n_components() corresponds to the dimension of the solution function in the PDE system, and correspondingly also to the number of equations your PDE system has. For example, the mixed Laplace system covered in step-20 has \(d+1\) components in \(d\) space dimensions: the scalar pressure and the \(d\) components of the velocity vector. Similarly, the elasticity equation covered in step-8 has \(d\) components in \(d\) space dimensions. In general, the number of components of a FESystem element is the accumulated number of components of all base elements times their multiplicities. A bit more on components is also given in the glossary entry on components.

While the concept of components is important from the viewpoint of a partial differential equation, the finite element side looks a bit different Since not only FESystem, but also vector-valued elements like FE_RaviartThomas, have several components. The concept needed here is a block. Each block encompasses the set of degrees of freedom associated with a single base element of an FESystem, where base elements with multiplicities count multiple times. These blocks are usually addressed using the information in DoFHandler::block_info(). The number of blocks of a FESystem object is simply the sum of all multiplicities of base elements and is given by n_blocks().

For example, the FESystem for the Taylor-Hood element for the three- dimensional Stokes problem can be built using the code

FE_Q<3> u(2);
FE_Q<3> p(1);
FESystem<3> sys1(u,3, p,1);

This example creates an FESystem sys1 with four components, three for the velocity components and one for the pressure, and also four blocks with the degrees of freedom of each of the velocity components and the pressure in a separate block each. The number of blocks is four since the first base element is repeated three times.

On the other hand, a Taylor-Hood element can also be constructed using

FESystem<3> sys2(U,1, p,1);

The FESystem sys2 created here has the same four components, but the degrees of freedom are distributed into only two blocks. The first block has all velocity degrees of freedom from U, while the second block contains the pressure degrees of freedom. Note that while U itself has 3 blocks, the FESystem sys2 does not attempt to split U into its base elements but considers it a block of its own. By blocking all velocities into one system first as in sys2, we achieve the same block structure that would be generated if instead of using a \(Q_2^3\) element for the velocities we had used vector-valued base elements, for instance like using a mixed discretization of Darcy's law using

FE_DGQ<3> p(1);
FESystem<3> sys3(u,1, p,1);

This example also produces a system with four components, but only two blocks.

In most cases, the composed element behaves as if it were a usual element. It just has more degrees of freedom than most of the "common" elements. However the underlying structure is visible in the restriction, prolongation and interface constraint matrices, which do not couple the degrees of freedom of the base elements. E.g. the continuity requirement is imposed for the shape functions of the subobjects separately; no requirement exist between shape functions of different subobjects, i.e. in the above example: on a hanging node, the respective value of the u velocity is only coupled to u at the vertices and the line on the larger cell next to this vertex, but there is no interaction with v and w of this or the other cell.

Internal information on numbering of degrees of freedom

The overall numbering of degrees of freedom is as follows: for each subobject (vertex, line, quad, or hex), the degrees of freedom are numbered such that we run over all subelements first, before turning for the next dof on this subobject or for the next subobject. For example, for an element of three components in one space dimension, the first two components being cubic lagrange elements and the third being a quadratic lagrange element, the ordering for the system s=(u,v,p) is:

That said, you should not rely on this numbering in your application as these internals might change in future. Rather use the functions system_to_component_index() and component_to_system_index().

For more information on the template parameter spacedim see the documentation of Triangulation.

Author
Wolfgang Bangerth, Guido Kanschat, 1999, 2002, 2003, 2006, Ralf Hartmann 2001.

Definition at line 44 of file fe.h.

Constructor & Destructor Documentation

◆ FESystem() [1/10]

template<int dim, int spacedim>
FESystem< dim, spacedim >::FESystem ( const FiniteElement< dim, spacedim > &  fe,
const unsigned int  n_elements 
)

Constructor. Take a finite element and the number of elements you want to group together using this class.

The object fe is not actually used for anything other than creating a copy that will then be owned by the current object. In other words, it is completely fine to call this constructor with a temporary object for the finite element, as in this code snippet:

Here, FE_Q<dim>(2) constructs an unnamed, temporary object that is passed to the FESystem constructor to create a finite element that consists of two components, both of which are quadratic FE_Q elements. The temporary is destroyed again at the end of the code that corresponds to this line, but this does not matter because FESystem creates its own copy of the FE_Q object.

This constructor (or its variants below) is used in essentially all tutorial programs that deal with vector valued problems. See step-8, step-20, step-22 and others for use cases. Also see the module on Handling vector valued problems.

Note
The material presented here is also discussed in video lecture 19, video lecture 20. (All video lectures are also available here.)
Parameters
[in]feThe finite element that will be used to represent the components of this composed element.
[in]n_elementsAn integer denoting how many copies of fe this element should consist of.

Definition at line 109 of file fe_system.cc.

◆ FESystem() [2/10]

template<int dim, int spacedim>
FESystem< dim, spacedim >::FESystem ( const FiniteElement< dim, spacedim > &  fe1,
const unsigned int  n1,
const FiniteElement< dim, spacedim > &  fe2,
const unsigned int  n2 
)

Constructor for mixed discretizations with two base elements.

See the other constructor above for an explanation of the general idea of composing elements.

Definition at line 128 of file fe_system.cc.

◆ FESystem() [3/10]

template<int dim, int spacedim>
FESystem< dim, spacedim >::FESystem ( const FiniteElement< dim, spacedim > &  fe1,
const unsigned int  n1,
const FiniteElement< dim, spacedim > &  fe2,
const unsigned int  n2,
const FiniteElement< dim, spacedim > &  fe3,
const unsigned int  n3 
)

Constructor for mixed discretizations with three base elements.

See the other constructor above for an explanation of the general idea of composing elements.

Definition at line 153 of file fe_system.cc.

◆ FESystem() [4/10]

template<int dim, int spacedim>
FESystem< dim, spacedim >::FESystem ( const FiniteElement< dim, spacedim > &  fe1,
const unsigned int  n1,
const FiniteElement< dim, spacedim > &  fe2,
const unsigned int  n2,
const FiniteElement< dim, spacedim > &  fe3,
const unsigned int  n3,
const FiniteElement< dim, spacedim > &  fe4,
const unsigned int  n4 
)

Constructor for mixed discretizations with four base elements.

See the first of the other constructors above for an explanation of the general idea of composing elements.

Definition at line 189 of file fe_system.cc.

◆ FESystem() [5/10]

template<int dim, int spacedim>
FESystem< dim, spacedim >::FESystem ( const FiniteElement< dim, spacedim > &  fe1,
const unsigned int  n1,
const FiniteElement< dim, spacedim > &  fe2,
const unsigned int  n2,
const FiniteElement< dim, spacedim > &  fe3,
const unsigned int  n3,
const FiniteElement< dim, spacedim > &  fe4,
const unsigned int  n4,
const FiniteElement< dim, spacedim > &  fe5,
const unsigned int  n5 
)

Constructor for mixed discretizations with five base elements.

See the first of the other constructors above for an explanation of the general idea of composing elements.

Definition at line 240 of file fe_system.cc.

◆ FESystem() [6/10]

template<int dim, int spacedim>
FESystem< dim, spacedim >::FESystem ( const std::vector< const FiniteElement< dim, spacedim > * > &  fes,
const std::vector< unsigned int > &  multiplicities 
)

Same as above but for any number of base elements. Pointers to the base elements and their multiplicities are passed as vectors to this constructor. The lengths of these vectors are assumed to be equal.

As above, the finite element objects pointed to by the first argument are not actually used other than to create copies internally. Consequently, you can delete these pointers immediately again after calling this constructor.

How to use this constructor

Using this constructor is a bit awkward at times because you need to pass two vectors in a place where it may not be straightforward to construct such a vector – for example, in the member initializer list of a class with an FESystem member variable. For example, if your main class looks like this:

template <int dim>
class MySimulator {
public:
MySimulator (const unsigned int polynomial_degree);
private:
};
template <int dim>
MySimulator<dim>::MySimulator (const unsigned int polynomial_degree)
:
fe (...) // what to pass here???
{}

Using the C++11 language standard (or later) you could do something like this to create an element with four base elements and multiplicities 1, 2, 3 and 4:

template <int dim>
MySimulator<dim>::MySimulator (const unsigned int polynomial_degree)
:
fe (std::vector<const FiniteElement<dim>*> { new FE_Q<dim>(1),
new FE_Q<dim>(2),
new FE_Q<dim>(3),
new FE_Q<dim>(4) },
std::vector<unsigned int> { 1, 2, 3, 4 })
{}

This creates two vectors in place and initializes them using the initializer list enclosed in braces { ... }.

This code has a problem: it creates four memory leaks because the first vector above is created with pointers to elements that are allocated with new but never destroyed.

The solution to the second of these problems is to create two static member functions that can create vectors. Here is an example:

template <int dim>
class MySimulator {
public:
MySimulator (const unsigned int polynomial_degree);
private:
static std::vector<const FiniteElement<dim>*>
create_fe_list (const unsigned int polynomial_degree);
static std::vector<unsigned int>
create_fe_multiplicities ();
};
template <int dim>
std::vector<const FiniteElement<dim>*>
MySimulator<dim>::create_fe_list (const unsigned int polynomial_degree)
{
std::vector<const FiniteElement<dim>*> fe_list;
fe_list.push_back (new FE_Q<dim>(1));
fe_list.push_back (new FE_Q<dim>(2));
fe_list.push_back (new FE_Q<dim>(3));
fe_list.push_back (new FE_Q<dim>(4));
return fe_list;
}
template <int dim>
std::vector<unsigned int>
MySimulator<dim>::create_fe_multiplicities ()
{
std::vector<unsigned int> multiplicities;
multiplicities.push_back (1);
multiplicities.push_back (2);
multiplicities.push_back (3);
multiplicities.push_back (4);
return multiplicities;
}
template <int dim>
MySimulator<dim>::MySimulator (const unsigned int polynomial_degree)
:
fe (create_fe_list (polynomial_degree),
create_fe_multiplicities ())
{}

The way this works is that we have two static member functions that create the necessary vectors to pass to the constructor of the member variable fe. They need to be static because they are called during the constructor of MySimulator at a time when the *this object isn't fully constructed and, consequently, regular member functions cannot be called yet.

The code above does not solve the problem with the memory leak yet, though: the create_fe_list() function creates a vector of pointers, but nothing destroys these. This is the solution:

template <int dim>
class MySimulator
{
public:
MySimulator (const unsigned int polynomial_degree);
private:
struct VectorElementDestroyer
{
const std::vector<const FiniteElement<dim>*> data;
VectorElementDestroyer(
const std::vector<const FiniteElement<dim>*> &pointers);
// destructor to delete the pointers
~VectorElementDestroyer ();
const std::vector<const FiniteElement<dim>*> & get_data () const;
};
static std::vector<const FiniteElement<dim>*>
create_fe_list (const unsigned int polynomial_degree);
static std::vector<unsigned int>
create_fe_multiplicities ();
};
template <int dim>
MySimulator<dim>::VectorElementDestroyer::
VectorElementDestroyer(
const std::vector<const FiniteElement<dim>*> &pointers)
:
data(pointers)
{}
template <int dim>
MySimulator<dim>::VectorElementDestroyer::
~VectorElementDestroyer ()
{
for (unsigned int i=0; i<data.size(); ++i)
delete data[i];
}
template <int dim>
const std::vector<const FiniteElement<dim>*> &
MySimulator<dim>::VectorElementDestroyer::
get_data () const
{
return data;
}
template <int dim>
MySimulator<dim>::MySimulator (const unsigned int polynomial_degree)
:
fe (VectorElementDestroyer(create_fe_list (polynomial_degree)).get_data(),
create_fe_multiplicities ())
{}

In other words, the vector we receive from the create_fe_list() is packed into a temporary object of type VectorElementDestroyer; we then get the vector from this temporary object immediately to pass it to the constructor of fe; and finally, the VectorElementDestroyer destructor is called at the end of the entire expression (after the constructor of fe has finished) and destroys the elements of the temporary vector. Voila: not short nor elegant, but it works!

Definition at line 293 of file fe_system.cc.

◆ FESystem() [7/10]

template<int dim, int spacedim = dim>
template<class... FEPairs, typename = typename enable_if_all< (std::is_same<typename std::decay<FEPairs>::type, std::pair<std::unique_ptr<FiniteElement<dim, spacedim>>, unsigned int>>::value || std::is_base_of<FiniteElement<dim, spacedim>, typename std::decay<FEPairs>::type>::value)...>::type>
FESystem< dim, spacedim >::FESystem ( FEPairs &&...  fe_pairs)

Constructor taking an arbitrary number of parameters of type std::pair<std::unique_ptr<FiniteElement<dim, spacedim>>, unsigned int>. In combination with FiniteElement::operator^, this allows to construct FESystem objects as follows:

FiniteElementType1<dim,spacedim> fe_1;
FiniteElementType1<dim,spacedim> fe_2;
FESystem<dim,spacedim> fe_system = ( fe_1^dim, fe_2^1 );

The FiniteElement objects are not actually used for anything other than creating a copy that will then be owned by the current object. In other words, it is completely fine to call this constructor with a temporary object for the finite element, as in this code snippet:

Here, FE_Q<dim>(2) constructs an unnamed, temporary object that is passed to the FESystem constructor to create a finite element that consists of two components, both of which are quadratic FE_Q elements. The temporary is destroyed again at the end of the code that corresponds to this line, but this does not matter because FESystem creates its own copy of the FE_Q object.

As a shortcut, this constructor also allows calling

instead of the more explicit

FESystem<dim> fe (FE_Q<dim>(2)^dim, FE_Q<dim>(1)^1);

In other words, if no multiplicity for an element is explicitly specified via the exponentiation operation, then it is assumed to be one (as one would have expected).

Warning
This feature is not available for Intel compilers prior to version 19.0. Defining this constructor leads to internal compiler errors for Intel compilers prior to 18.0.

◆ FESystem() [8/10]

template<int dim, int spacedim = dim>
FESystem< dim, spacedim >::FESystem ( const std::initializer_list< std::pair< std::unique_ptr< FiniteElement< dim, spacedim >>, unsigned int >> &  fe_systems)

Same as above allowing the following syntax:

FiniteElementType1<dim,spacedim> fe_1;
FiniteElementType1<dim,spacedim> fe_2;
FESystem<dim,spacedim> fe_system = { fe_1^dim, fe_2^1 };
Warning
This feature is not available for Intel compilers prior to version 19.0. The constructor is just not selected for overload resolution.

◆ FESystem() [9/10]

template<int dim, int spacedim = dim>
FESystem< dim, spacedim >::FESystem ( const FESystem< dim, spacedim > &  )
delete

Copy constructor. This constructor is deleted, i.e., copying FESystem objects is not allowed.

◆ FESystem() [10/10]

template<int dim, int spacedim = dim>
FESystem< dim, spacedim >::FESystem ( FESystem< dim, spacedim > &&  other_fe_system)
inlinenoexcept

Move constructor.

Definition at line 532 of file fe_system.h.

◆ ~FESystem()

template<int dim, int spacedim = dim>
virtual FESystem< dim, spacedim >::~FESystem ( )
overridevirtualdefault

Destructor.

Member Function Documentation

◆ get_name()

template<int dim, int spacedim>
std::string FESystem< dim, spacedim >::get_name
overridevirtual

Return a string that uniquely identifies a finite element. This element returns a string that is composed of the strings name1...nameN returned by the basis elements. From these, we create a sequence FESystem<dim>[name1^m1-name2^m2-...-nameN^mN], where mi are the multiplicities of the basis elements. If a multiplicity is equal to one, then the superscript is omitted.

Definition at line 311 of file fe_system.cc.

◆ clone()

template<int dim, int spacedim>
std::unique_ptr< FiniteElement< dim, spacedim > > FESystem< dim, spacedim >::clone
overridevirtual

Definition at line 340 of file fe_system.cc.

◆ requires_update_flags()

template<int dim, int spacedim>
UpdateFlags FESystem< dim, spacedim >::requires_update_flags ( const UpdateFlags  update_flags) const
overridevirtual

Definition at line 896 of file fe_system.cc.

◆ get_sub_fe()

template<int dim, int spacedim>
const FiniteElement< dim, spacedim > & FESystem< dim, spacedim >::get_sub_fe ( const unsigned int  first_component,
const unsigned int  n_selected_components 
) const
overridevirtual

Return a reference to a contained finite element that matches the components selected by the given ComponentMask mask.

For an arbitrarily nested FESystem, this function returns the inner-most FiniteElement that matches the given mask. The method fails if the mask does not exactly match one of the contained finite elements. It is most useful if the current object is an FESystem, as the return value can only be this in all other cases.

Note that the returned object can be an FESystem if the mask matches it but not any of the contained objects.

Let us illustrate the function with the an FESystem fe with 7 components:

FESystem<2> fe_velocity(FE_Q<2>(2), 2);
FE_Q<2> fe_pressure(1);
FE_DGP<2> fe_dg(0);
FE_BDM<2> fe_nonprim(1);
FESystem<2> fe(fe_velocity, 1, fe_pressure, 1, fe_dg, 2, fe_nonprim, 1);

The following table lists all possible component masks you can use:

ComponentMask Result Description
[true,true,true,true,true,true,true] FESystem<2>[FESystem<2>[FE_Q<2>(2)^2]-FE_Q<2>(1)-FE_DGP<2>(0)^2-FE_BDM<2>(1)] fe itself, the whole FESystem
[true,true,false,false,false,false,false] FESystem<2>[FE_Q<2>(2)^2] just the fe_velocity
[true,false,false,false,false,false,false] FE_Q<2>(2) The first component in fe_velocity
[false,true,false,false,false,false,false] FE_Q<2>(2) The second component in fe_velocity
[false,false,true,false,false,false,false] FE_Q<2>(1) fe_pressure
[false,false,false,true,false,false,false] FE_DGP<2>(0) first copy of fe_dg
[false,false,false,false,true,false,false] FE_DGP<2>(0) second copy of fe_dg
[false,false,false,false,false,true,true] FE_BDM<2>(1) both components of fe_nonprim

Definition at line 357 of file fe_system.cc.

◆ shape_value()

template<int dim, int spacedim>
double FESystem< dim, spacedim >::shape_value ( const unsigned int  i,
const Point< dim > &  p 
) const
overridevirtual

Return the value of the ith shape function at the point p. p is a point on the reference element. Since this finite element is always vector-valued, we return the value of the only non-zero component of the vector value of this shape function. If the shape function has more than one non-zero component (which we refer to with the term non-primitive), then throw an exception of type ExcShapeFunctionNotPrimitive.

An ExcUnitShapeValuesDoNotExist is thrown if the shape values of the FiniteElement (corresponding to the ith shape function) depend on the shape of the cell in real space.

Definition at line 386 of file fe_system.cc.

◆ shape_value_component()

template<int dim, int spacedim>
double FESystem< dim, spacedim >::shape_value_component ( const unsigned int  i,
const Point< dim > &  p,
const unsigned int  component 
) const
overridevirtual

Return the value of the componentth vector component of the ith shape function at the point p. See the FiniteElement base class for more information about the semantics of this function.

Since this element is vector valued in general, it relays the computation of these values to the base elements.

Definition at line 402 of file fe_system.cc.

◆ shape_grad()

template<int dim, int spacedim>
Tensor< 1, dim > FESystem< dim, spacedim >::shape_grad ( const unsigned int  i,
const Point< dim > &  p 
) const
overridevirtual

Return the gradient of the ith shape function at the point p. p is a point on the reference element, and likewise the gradient is the gradient on the unit cell with respect to unit cell coordinates. Since this finite element is always vector-valued, we return the value of the only non-zero component of the vector value of this shape function. If the shape function has more than one non-zero component (which we refer to with the term non-primitive), then throw an exception of type ExcShapeFunctionNotPrimitive.

An ExcUnitShapeValuesDoNotExist is thrown if the shape values of the FiniteElement (corresponding to the ith shape function) depend on the shape of the cell in real space.

Definition at line 438 of file fe_system.cc.

◆ shape_grad_component()

template<int dim, int spacedim>
Tensor< 1, dim > FESystem< dim, spacedim >::shape_grad_component ( const unsigned int  i,
const Point< dim > &  p,
const unsigned int  component 
) const
overridevirtual

Return the gradient of the componentth vector component of the ith shape function at the point p. See the FiniteElement base class for more information about the semantics of this function.

Since this element is vector valued in general, it relays the computation of these values to the base elements.

Definition at line 454 of file fe_system.cc.

◆ shape_grad_grad()

template<int dim, int spacedim>
Tensor< 2, dim > FESystem< dim, spacedim >::shape_grad_grad ( const unsigned int  i,
const Point< dim > &  p 
) const
overridevirtual

Return the tensor of second derivatives of the ith shape function at point p on the unit cell. The derivatives are derivatives on the unit cell with respect to unit cell coordinates. Since this finite element is always vector-valued, we return the value of the only non-zero component of the vector value of this shape function. If the shape function has more than one non-zero component (which we refer to with the term non- primitive), then throw an exception of type ExcShapeFunctionNotPrimitive.

An ExcUnitShapeValuesDoNotExist is thrown if the shape values of the FiniteElement (corresponding to the ith shape function) depend on the shape of the cell in real space.

Definition at line 483 of file fe_system.cc.

◆ shape_grad_grad_component()

template<int dim, int spacedim>
Tensor< 2, dim > FESystem< dim, spacedim >::shape_grad_grad_component ( const unsigned int  i,
const Point< dim > &  p,
const unsigned int  component 
) const
overridevirtual

Return the second derivatives of the componentth vector component of the ith shape function at the point p. See the FiniteElement base class for more information about the semantics of this function.

Since this element is vector valued in general, it relays the computation of these values to the base elements.

Definition at line 499 of file fe_system.cc.

◆ shape_3rd_derivative()

template<int dim, int spacedim>
Tensor< 3, dim > FESystem< dim, spacedim >::shape_3rd_derivative ( const unsigned int  i,
const Point< dim > &  p 
) const
overridevirtual

Return the tensor of third derivatives of the ith shape function at point p on the unit cell. The derivatives are derivatives on the unit cell with respect to unit cell coordinates. Since this finite element is always vector-valued, we return the value of the only non-zero component of the vector value of this shape function. If the shape function has more than one non-zero component (which we refer to with the term non- primitive), then throw an exception of type ExcShapeFunctionNotPrimitive.

An ExcUnitShapeValuesDoNotExist is thrown if the shape values of the FiniteElement (corresponding to the ith shape function) depend on the shape of the cell in real space.

Definition at line 528 of file fe_system.cc.

◆ shape_3rd_derivative_component()

template<int dim, int spacedim>
Tensor< 3, dim > FESystem< dim, spacedim >::shape_3rd_derivative_component ( const unsigned int  i,
const Point< dim > &  p,
const unsigned int  component 
) const
overridevirtual

Return the third derivatives of the componentth vector component of the ith shape function at the point p. See the FiniteElement base class for more information about the semantics of this function.

Since this element is vector valued in general, it relays the computation of these values to the base elements.

Definition at line 544 of file fe_system.cc.

◆ shape_4th_derivative()

template<int dim, int spacedim>
Tensor< 4, dim > FESystem< dim, spacedim >::shape_4th_derivative ( const unsigned int  i,
const Point< dim > &  p 
) const
overridevirtual

Return the tensor of fourth derivatives of the ith shape function at point p on the unit cell. The derivatives are derivatives on the unit cell with respect to unit cell coordinates. Since this finite element is always vector-valued, we return the value of the only non-zero component of the vector value of this shape function. If the shape function has more than one non-zero component (which we refer to with the term non- primitive), then throw an exception of type ExcShapeFunctionNotPrimitive.

An ExcUnitShapeValuesDoNotExist is thrown if the shape values of the FiniteElement (corresponding to the ith shape function) depend on the shape of the cell in real space.

Definition at line 573 of file fe_system.cc.

◆ shape_4th_derivative_component()

template<int dim, int spacedim>
Tensor< 4, dim > FESystem< dim, spacedim >::shape_4th_derivative_component ( const unsigned int  i,
const Point< dim > &  p,
const unsigned int  component 
) const
overridevirtual

Return the fourth derivatives of the componentth vector component of the ith shape function at the point p. See the FiniteElement base class for more information about the semantics of this function.

Since this element is vector valued in general, it relays the computation of these values to the base elements.

Definition at line 589 of file fe_system.cc.

◆ get_interpolation_matrix()

template<int dim, int spacedim>
void FESystem< dim, spacedim >::get_interpolation_matrix ( const FiniteElement< dim, spacedim > &  source,
FullMatrix< double > &  matrix 
) const
overridevirtual

Return the matrix interpolating from the given finite element to the present one. The size of the matrix is then dofs_per_cell times source.dofs_per_cell.

These matrices are available if source and destination element are both FESystem elements, have the same number of base elements with same element multiplicity, and if these base elements also implement their get_interpolation_matrix functions. Otherwise, an exception of type FiniteElement<dim,spacedim>::ExcInterpolationNotImplemented is thrown.

Definition at line 618 of file fe_system.cc.

◆ base_element()

template<int dim, int spacedim>
const FiniteElement< dim, spacedim > & FESystem< dim, spacedim >::base_element ( const unsigned int  index) const
overridevirtual

Access to a composing element. The index needs to be smaller than the number of base elements. Note that the number of base elements may in turn be smaller than the number of components of the system element, if the multiplicities are greater than one.

Definition at line 2286 of file fe_system.cc.

◆ has_support_on_face()

template<int dim, int spacedim>
bool FESystem< dim, spacedim >::has_support_on_face ( const unsigned int  shape_index,
const unsigned int  face_index 
) const
overridevirtual

This function returns true, if the shape function shape_index has non-zero function values somewhere on the face face_index.

Definition at line 2296 of file fe_system.cc.

◆ get_restriction_matrix()

template<int dim, int spacedim>
const FullMatrix< double > & FESystem< dim, spacedim >::get_restriction_matrix ( const unsigned int  child,
const RefinementCase< dim > &  refinement_case = RefinementCase<dim>::isotropic_refinement 
) const
overridevirtual

Projection from a fine grid space onto a coarse grid space. Overrides the respective method in FiniteElement, implementing lazy evaluation (initialize when requested).

If this projection operator is associated with a matrix P, then the restriction of this matrix P_i to a single child cell is returned here.

The matrix P is the concatenation or the sum of the cell matrices P_i, depending on the #restriction_is_additive_flags. This distinguishes interpolation (concatenation) and projection with respect to scalar products (summation).

Row and column indices are related to coarse grid and fine grid spaces, respectively, consistent with the definition of the associated operator.

If projection matrices are not implemented in the derived finite element class, this function aborts with an exception of type FiniteElement::ExcProjectionVoid. You can check whether this would happen by first calling the restriction_is_implemented() or the isotropic_restriction_is_implemented() function.

Definition at line 694 of file fe_system.cc.

◆ get_prolongation_matrix()

template<int dim, int spacedim>
const FullMatrix< double > & FESystem< dim, spacedim >::get_prolongation_matrix ( const unsigned int  child,
const RefinementCase< dim > &  refinement_case = RefinementCase<dim>::isotropic_refinement 
) const
overridevirtual

Embedding matrix between grids. Overrides the respective method in FiniteElement, implementing lazy evaluation (initialize when queried).

The identity operator from a coarse grid space into a fine grid space is associated with a matrix P. The restriction of this matrix P_i to a single child cell is returned here.

The matrix P is the concatenation, not the sum of the cell matrices P_i. That is, if the same non-zero entry j,k exists in two different child matrices P_i, the value should be the same in both matrices and it is copied into the matrix P only once.

Row and column indices are related to fine grid and coarse grid spaces, respectively, consistent with the definition of the associated operator.

These matrices are used by routines assembling the prolongation matrix for multi-level methods. Upon assembling the transfer matrix between cells using this matrix array, zero elements in the prolongation matrix are discarded and will not fill up the transfer matrix.

If prolongation matrices are not implemented in one of the base finite element classes, this function aborts with an exception of type FiniteElement::ExcEmbeddingVoid. You can check whether this would happen by first calling the prolongation_is_implemented() or the isotropic_prolongation_is_implemented() function.

Definition at line 785 of file fe_system.cc.

◆ face_to_cell_index()

template<int dim, int spacedim>
unsigned int FESystem< dim, spacedim >::face_to_cell_index ( const unsigned int  face_dof_index,
const unsigned int  face,
const bool  face_orientation = true,
const bool  face_flip = false,
const bool  face_rotation = false 
) const
overridevirtual

Given an index in the natural ordering of indices on a face, return the index of the same degree of freedom on the cell.

To explain the concept, consider the case where we would like to know whether a degree of freedom on a face, for example as part of an FESystem element, is primitive. Unfortunately, the is_primitive() function in the FiniteElement class takes a cell index, so we would need to find the cell index of the shape function that corresponds to the present face index. This function does that.

Code implementing this would then look like this:

for (i=0; i<dofs_per_face; ++i)
if (fe.is_primitive(fe.face_to_cell_index(i, some_face_no)))
... do whatever

The function takes additional arguments that account for the fact that actual faces can be in their standard ordering with respect to the cell under consideration, or can be flipped, oriented, etc.

Parameters
face_dof_indexThe index of the degree of freedom on a face. This index must be between zero and dofs_per_face.
faceThe number of the face this degree of freedom lives on. This number must be between zero and GeometryInfo::faces_per_cell.
face_orientationOne part of the description of the orientation of the face. See GlossFaceOrientation.
face_flipOne part of the description of the orientation of the face. See GlossFaceOrientation.
face_rotationOne part of the description of the orientation of the face. See GlossFaceOrientation.
Returns
The index of this degree of freedom within the set of degrees of freedom on the entire cell. The returned value will be between zero and dofs_per_cell.

Definition at line 851 of file fe_system.cc.

◆ unit_support_point()

template<int dim, int spacedim>
Point< dim > FESystem< dim, spacedim >::unit_support_point ( const unsigned int  index) const
overridevirtual

Implementation of the respective function in the base class.

Definition at line 2309 of file fe_system.cc.

◆ unit_face_support_point()

template<int dim, int spacedim>
Point< dim - 1 > FESystem< dim, spacedim >::unit_face_support_point ( const unsigned int  index) const
overridevirtual

Implementation of the respective function in the base class.

Definition at line 2330 of file fe_system.cc.

◆ get_constant_modes()

template<int dim, int spacedim>
std::pair< Table< 2, bool >, std::vector< unsigned int > > FESystem< dim, spacedim >::get_constant_modes
overridevirtual

Return a list of constant modes of the element. The returns table has as many rows as there are components in the element and dofs_per_cell columns. To each component of the finite element, the row in the returned table contains a basis representation of the constant function 1 on the element. Concatenates the constant modes of each base element.

Definition at line 2352 of file fe_system.cc.

◆ hp_constraints_are_implemented()

template<int dim, int spacedim>
bool FESystem< dim, spacedim >::hp_constraints_are_implemented
overridevirtual

Return whether this element implements its hanging node constraints in the new way, which has to be used to make elements "hp compatible".

This function returns true if and only if all its base elements return true for this function.

Definition at line 1858 of file fe_system.cc.

◆ get_face_interpolation_matrix()

template<int dim, int spacedim>
void FESystem< dim, spacedim >::get_face_interpolation_matrix ( const FiniteElement< dim, spacedim > &  source,
FullMatrix< double > &  matrix 
) const
overridevirtual

Return the matrix interpolating from a face of one element to the face of the neighboring element. The size of the matrix is then source.dofs_per_face times this->dofs_per_face.

Base elements of this element will have to implement this function. They may only provide interpolation matrices for certain source finite elements, for example those from the same family. If they don't implement interpolation from a given element, then they must throw an exception of type FiniteElement<dim,spacedim>::ExcInterpolationNotImplemented, which will get propagated out from this element.

Definition at line 1871 of file fe_system.cc.

◆ get_subface_interpolation_matrix()

template<int dim, int spacedim>
void FESystem< dim, spacedim >::get_subface_interpolation_matrix ( const FiniteElement< dim, spacedim > &  source,
const unsigned int  subface,
FullMatrix< double > &  matrix 
) const
overridevirtual

Return the matrix interpolating from a face of one element to the subface of the neighboring element. The size of the matrix is then source.dofs_per_face times this->dofs_per_face.

Base elements of this element will have to implement this function. They may only provide interpolation matrices for certain source finite elements, for example those from the same family. If they don't implement interpolation from a given element, then they must throw an exception of type FiniteElement<dim,spacedim>::ExcInterpolationNotImplemented, which will get propagated out from this element.

Definition at line 1980 of file fe_system.cc.

◆ hp_vertex_dof_identities()

template<int dim, int spacedim>
std::vector< std::pair< unsigned int, unsigned int > > FESystem< dim, spacedim >::hp_vertex_dof_identities ( const FiniteElement< dim, spacedim > &  fe_other) const
overridevirtual

If, on a vertex, several finite elements are active, the hp code first assigns the degrees of freedom of each of these FEs different global indices. It then calls this function to find out which of them should get identical values, and consequently can receive the same global DoF index. This function therefore returns a list of identities between DoFs of the present finite element object with the DoFs of fe_other, which is a reference to a finite element object representing one of the other finite elements active on this particular vertex. The function computes which of the degrees of freedom of the two finite element objects are equivalent, both numbered between zero and the corresponding value of dofs_per_vertex of the two finite elements. The first index of each pair denotes one of the vertex dofs of the present element, whereas the second is the corresponding index of the other finite element.

Definition at line 2208 of file fe_system.cc.

◆ hp_line_dof_identities()

template<int dim, int spacedim>
std::vector< std::pair< unsigned int, unsigned int > > FESystem< dim, spacedim >::hp_line_dof_identities ( const FiniteElement< dim, spacedim > &  fe_other) const
overridevirtual

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

Definition at line 2216 of file fe_system.cc.

◆ hp_quad_dof_identities()

template<int dim, int spacedim>
std::vector< std::pair< unsigned int, unsigned int > > FESystem< dim, spacedim >::hp_quad_dof_identities ( const FiniteElement< dim, spacedim > &  fe_other) const
overridevirtual

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

Definition at line 2226 of file fe_system.cc.

◆ compare_for_domination()

template<int dim, int spacedim>
FiniteElementDomination::Domination FESystem< dim, spacedim >::compare_for_domination ( const FiniteElement< dim, spacedim > &  fe_other,
const unsigned int  codim = 0 
) const
finaloverridevirtual

Return whether this element dominates another one given as argument fe_other, whether it is the other way around, whether neither dominates, or if either could dominate. The codim parameter describes the codimension of the investigated subspace and specifies that it is subject to this comparison. For example, if codim==0 then this function compares which element dominates at the cell level. If codim==1, then the elements are compared at faces, i.e., the comparison happens between the function spaces of the two finite elements as restricted to a face. Larger values of codim work correspondingly.

For a definition of domination, see FiniteElementDomination::Domination and in particular the hp paper.

Definition at line 2236 of file fe_system.cc.

◆ convert_generalized_support_point_values_to_dof_values()

template<int dim, int spacedim>
void FESystem< dim, spacedim >::convert_generalized_support_point_values_to_dof_values ( const std::vector< Vector< double >> &  support_point_values,
std::vector< double > &  dof_values 
) const
overridevirtual

Implementation of the FiniteElement::convert_generalized_support_point_values_to_dof_values() function.

This function simply calls FiniteElement::convert_generalized_support_point_values_to_dof_values of the base elements and re-assembles everything into the output argument. If a base element is non-interpolatory the corresponding dof values are filled with "signaling" NaNs instead.

The function fails if none of the base elements of the FESystem are interpolatory.

Definition at line 2408 of file fe_system.cc.

◆ memory_consumption()

template<int dim, int spacedim>
std::size_t FESystem< dim, spacedim >::memory_consumption
overridevirtual

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

This function is made virtual, since finite element objects are usually accessed through pointers to their base class, rather than the class itself.

Definition at line 2515 of file fe_system.cc.

◆ get_data()

template<int dim, int spacedim>
std::unique_ptr< typename FiniteElement< dim, spacedim >::InternalDataBase > FESystem< dim, spacedim >::get_data ( const UpdateFlags  update_flags,
const Mapping< dim, spacedim > &  mapping,
const Quadrature< dim > &  quadrature,
::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &  output_data 
) const
overrideprotectedvirtual

Definition at line 910 of file fe_system.cc.

◆ get_face_data()

template<int dim, int spacedim>
std::unique_ptr< typename FiniteElement< dim, spacedim >::InternalDataBase > FESystem< dim, spacedim >::get_face_data ( const UpdateFlags  update_flags,
const Mapping< dim, spacedim > &  mapping,
const Quadrature< dim - 1 > &  quadrature,
::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &  output_data 
) const
overrideprotectedvirtual

Definition at line 973 of file fe_system.cc.

◆ get_subface_data()

template<int dim, int spacedim>
std::unique_ptr< typename FiniteElement< dim, spacedim >::InternalDataBase > FESystem< dim, spacedim >::get_subface_data ( const UpdateFlags  update_flags,
const Mapping< dim, spacedim > &  mapping,
const Quadrature< dim - 1 > &  quadrature,
::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &  output_data 
) const
overrideprotectedvirtual

Definition at line 1036 of file fe_system.cc.

◆ fill_fe_values()

template<int dim, int spacedim>
void FESystem< dim, spacedim >::fill_fe_values ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell,
const CellSimilarity::Similarity  cell_similarity,
const Quadrature< dim > &  quadrature,
const Mapping< dim, spacedim > &  mapping,
const typename Mapping< dim, spacedim >::InternalDataBase &  mapping_internal,
const ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &  mapping_data,
const typename FiniteElement< dim, spacedim >::InternalDataBase &  fe_internal,
::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &  output_data 
) const
overrideprotectedvirtual

Definition at line 1097 of file fe_system.cc.

◆ fill_fe_face_values()

template<int dim, int spacedim>
void FESystem< dim, spacedim >::fill_fe_face_values ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell,
const unsigned int  face_no,
const Quadrature< dim - 1 > &  quadrature,
const Mapping< dim, spacedim > &  mapping,
const typename Mapping< dim, spacedim >::InternalDataBase &  mapping_internal,
const ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &  mapping_data,
const typename FiniteElement< dim, spacedim >::InternalDataBase &  fe_internal,
::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &  output_data 
) const
overrideprotectedvirtual

Definition at line 1127 of file fe_system.cc.

◆ fill_fe_subface_values()

template<int dim, int spacedim>
void FESystem< dim, spacedim >::fill_fe_subface_values ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell,
const unsigned int  face_no,
const unsigned int  sub_no,
const Quadrature< dim - 1 > &  quadrature,
const Mapping< dim, spacedim > &  mapping,
const typename Mapping< dim, spacedim >::InternalDataBase &  mapping_internal,
const ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &  mapping_data,
const typename FiniteElement< dim, spacedim >::InternalDataBase &  fe_internal,
::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &  output_data 
) const
overrideprotectedvirtual

Definition at line 1157 of file fe_system.cc.

◆ compute_fill()

template<int dim, int spacedim>
template<int dim_1>
void FESystem< dim, spacedim >::compute_fill ( const Mapping< dim, spacedim > &  mapping,
const typename Triangulation< dim, spacedim >::cell_iterator &  cell,
const unsigned int  face_no,
const unsigned int  sub_no,
const Quadrature< dim_1 > &  quadrature,
const CellSimilarity::Similarity  cell_similarity,
const typename Mapping< dim, spacedim >::InternalDataBase &  mapping_internal,
const typename FiniteElement< dim, spacedim >::InternalDataBase &  fe_data,
const internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &  mapping_data,
internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &  output_data 
) const
protected

Do the work for the three fill_fe*_values functions.

Calls (among other things) fill_fe_([sub]face)_values of the base elements. Calls fill_fe_values if face_no==invalid_face_no and sub_no==invalid_face_no; calls fill_fe_face_values if face_no==invalid_face_no and sub_no!=invalid_face_no; and calls fill_fe_subface_values if face_no!=invalid_face_no and sub_no!=invalid_face_no.

Definition at line 1189 of file fe_system.cc.

◆ initialize()

template<int dim, int spacedim>
void FESystem< dim, spacedim >::initialize ( const std::vector< const FiniteElement< dim, spacedim > * > &  fes,
const std::vector< unsigned int > &  multiplicities 
)
private

This function is simply singled out of the constructors since there are several of them. It sets up the index table for the system as well as restriction and prolongation matrices.

Definition at line 1595 of file fe_system.cc.

◆ build_interface_constraints()

template<int dim, int spacedim>
void FESystem< dim, spacedim >::build_interface_constraints
private

Used by initialize.

Definition at line 1391 of file fe_system.cc.

◆ hp_object_dof_identities()

template<int dim, int spacedim>
template<int structdim>
std::vector< std::pair< unsigned int, unsigned int > > FESystem< dim, spacedim >::hp_object_dof_identities ( const FiniteElement< dim, spacedim > &  fe_other) const
private

A function that computes the hp_vertex_dof_identities(), hp_line_dof_identities(), or hp_quad_dof_identities(), depending on the value of the template parameter.

Definition at line 2097 of file fe_system.cc.

Friends And Related Function Documentation

◆ FE_Enriched< dim, spacedim >

template<int dim, int spacedim = dim>
friend class FE_Enriched< dim, spacedim >
friend

Definition at line 1242 of file fe_system.h.

Member Data Documentation

◆ invalid_face_number

template<int dim, int spacedim>
const unsigned int FESystem< dim, spacedim >::invalid_face_number = numbers::invalid_unsigned_int
staticprivate

Value to indicate that a given face or subface number is invalid.

Definition at line 1113 of file fe_system.h.

◆ base_elements

template<int dim, int spacedim = dim>
std::vector<std::pair<std::unique_ptr<const FiniteElement<dim, spacedim> >, unsigned int> > FESystem< dim, spacedim >::base_elements
private

Pointers to underlying finite element objects.

This object contains a pointer to each contributing element of a mixed discretization and its multiplicity. It is created by the constructor and constant afterwards.

Definition at line 1124 of file fe_system.h.

◆ generalized_support_points_index_table

template<int dim, int spacedim = dim>
std::vector<std::vector<std::size_t> > FESystem< dim, spacedim >::generalized_support_points_index_table
private

An index table that maps generalized support points of a base element to the vector of generalized support points of the FE System. It holds true that

generalized_support_points[n] ==
base_elements[i].generalized_support_points[j];

for each base element (indexed by i) and each g. s. point of the base element (index by j).

Definition at line 1138 of file fe_system.h.

◆ mutex

template<int dim, int spacedim = dim>
std::mutex FESystem< dim, spacedim >::mutex
mutableprivate

Mutex for protecting initialization of restriction and embedding matrix.

Definition at line 1240 of file fe_system.h.


The documentation for this class was generated from the following files:
FE_DGQ
Definition: fe_dgq.h:112
FE_Q
Definition: fe_q.h:554
FE_RaviartThomas
Definition: fe_raviart_thomas.h:130
LAPACKSupport::U
static const char U
Definition: lapack_support.h:167
FESystem::base_elements
std::vector< std::pair< std::unique_ptr< const FiniteElement< dim, spacedim > >, unsigned int > > base_elements
Definition: fe_system.h:1124
FESystem::get_data
virtual std::unique_ptr< typename FiniteElement< dim, spacedim >::InternalDataBase > get_data(const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
Definition: fe_system.cc:910
FESystem::generalized_support_points_index_table
std::vector< std::vector< std::size_t > > generalized_support_points_index_table
Definition: fe_system.h:1138
FE_BDM
Definition: fe_bdm.h:59
FiniteElement
Definition: fe.h:648
FE_DGP
Definition: fe_dgp.h:311
FESystem
Definition: fe.h:44