Reference documentation for deal.II version 9.0.0
Public Member Functions | Static Private Member Functions | Private Attributes | List of all members
FE_Q_Bubbles< dim, spacedim > Class Template Reference

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

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

Public Member Functions

 FE_Q_Bubbles (const unsigned int p)
 
 FE_Q_Bubbles (const Quadrature< 1 > &points)
 
virtual std::string get_name () const
 
virtual void convert_generalized_support_point_values_to_dof_values (const std::vector< Vector< double > > &support_point_values, std::vector< double > &nodal_values) const
 
virtual void get_interpolation_matrix (const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
 
virtual const FullMatrix< double > & get_prolongation_matrix (const unsigned int child, const RefinementCase< dim > &refinement_case) const
 
virtual const FullMatrix< double > & get_restriction_matrix (const unsigned int child, const RefinementCase< dim > &refinement_case) const
 
virtual bool has_support_on_face (const unsigned int shape_index, const unsigned int face_index) const
 
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone () const
 
- Public Member Functions inherited from FE_Q_Base< TensorProductPolynomialsBubbles< dim >, dim, spacedim >
 FE_Q_Base (const TensorProductPolynomialsBubbles< dim > &poly_space, const FiniteElementData< dim > &fe_data, const std::vector< bool > &restriction_is_additive_flags)
 
virtual void get_face_interpolation_matrix (const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
 
virtual void get_subface_interpolation_matrix (const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix) const
 
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
 
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes () const
 
virtual bool hp_constraints_are_implemented () const
 
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities (const FiniteElement< dim, spacedim > &fe_other) const
 
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities (const FiniteElement< dim, spacedim > &fe_other) const
 
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_quad_dof_identities (const FiniteElement< dim, spacedim > &fe_other) const
 
virtual FiniteElementDomination::Domination compare_for_face_domination (const FiniteElement< dim, spacedim > &fe_other) const
 
- Public Member Functions inherited from FE_Poly< TensorProductPolynomialsBubbles< dim >, dim, spacedim >
 FE_Poly (const TensorProductPolynomialsBubbles< dim > &poly_space, const FiniteElementData< dim > &fe_data, const std::vector< bool > &restriction_is_additive_flags, const std::vector< ComponentMask > &nonzero_components)
 
unsigned int get_degree () const
 
std::vector< unsigned int > get_poly_space_numbering () const
 
std::vector< unsigned int > get_poly_space_numbering_inverse () const
 
virtual double shape_value (const unsigned int i, const Point< dim > &p) const
 
virtual double shape_value_component (const unsigned int i, const Point< dim > &p, const unsigned int component) const
 
virtual Tensor< 1, dim > shape_grad (const unsigned int i, const Point< dim > &p) const
 
virtual Tensor< 1, dim > shape_grad_component (const unsigned int i, const Point< dim > &p, const unsigned int component) const
 
virtual Tensor< 2, dim > shape_grad_grad (const unsigned int i, const Point< dim > &p) const
 
virtual Tensor< 2, dim > shape_grad_grad_component (const unsigned int i, const Point< dim > &p, const unsigned int component) const
 
virtual Tensor< 3, dim > shape_3rd_derivative (const unsigned int i, const Point< dim > &p) const
 
virtual Tensor< 3, dim > shape_3rd_derivative_component (const unsigned int i, const Point< dim > &p, const unsigned int component) const
 
virtual Tensor< 4, dim > shape_4th_derivative (const unsigned int i, const Point< dim > &p) const
 
virtual Tensor< 4, dim > shape_4th_derivative_component (const unsigned int i, const Point< dim > &p, const unsigned int component) const
 
- Public Member Functions inherited from FiniteElement< dim, spacedim >
 FiniteElement (const FiniteElementData< dim > &fe_data, const std::vector< bool > &restriction_is_additive_flags, const std::vector< ComponentMask > &nonzero_components)
 
 FiniteElement (FiniteElement< dim, spacedim > &&)=default
 
 FiniteElement (const FiniteElement< dim, spacedim > &)=default
 
virtual ~FiniteElement ()=default
 
std::pair< std::unique_ptr< FiniteElement< dim, spacedim > >, unsigned int > operator^ (const unsigned int multiplicity) const
 
const FiniteElement< dim, spacedim > & operator[] (const unsigned int fe_index) const
 
virtual bool operator== (const FiniteElement< dim, spacedim > &fe) const
 
bool operator!= (const FiniteElement< dim, spacedim > &) const
 
virtual std::size_t memory_consumption () const
 
bool prolongation_is_implemented () const
 
bool isotropic_prolongation_is_implemented () const
 
bool restriction_is_implemented () const
 
bool isotropic_restriction_is_implemented () const
 
bool restriction_is_additive (const unsigned int index) const
 
const FullMatrix< double > & constraints (const ::internal::SubfaceCase< dim > &subface_case=::internal::SubfaceCase< dim >::case_isotropic) const
 
bool constraints_are_implemented (const ::internal::SubfaceCase< dim > &subface_case=::internal::SubfaceCase< dim >::case_isotropic) const
 
std::pair< unsigned int, unsigned int > system_to_component_index (const unsigned int index) const
 
unsigned int component_to_system_index (const unsigned int component, const unsigned int index) const
 
std::pair< unsigned int, unsigned int > face_system_to_component_index (const unsigned int index) const
 
unsigned int adjust_quad_dof_index_for_face_orientation (const unsigned int index, const bool face_orientation, const bool face_flip, const bool face_rotation) const
 
unsigned int adjust_line_dof_index_for_line_orientation (const unsigned int index, const bool line_orientation) const
 
const ComponentMaskget_nonzero_components (const unsigned int i) const
 
unsigned int n_nonzero_components (const unsigned int i) const
 
bool is_primitive () const
 
bool is_primitive (const unsigned int i) const
 
unsigned int n_base_elements () const
 
virtual const FiniteElement< dim, spacedim > & base_element (const unsigned int index) const
 
unsigned int element_multiplicity (const unsigned int index) const
 
const FiniteElement< dim, spacedim > & get_sub_fe (const ComponentMask &mask) const
 
virtual const FiniteElement< dim, spacedim > & get_sub_fe (const unsigned int first_component, const unsigned int n_selected_components) const
 
std::pair< std::pair< unsigned int, unsigned int >, unsigned int > system_to_base_index (const unsigned int index) const
 
std::pair< std::pair< unsigned int, unsigned int >, unsigned int > face_system_to_base_index (const unsigned int index) const
 
types::global_dof_index first_block_of_base (const unsigned int b) const
 
std::pair< unsigned int, unsigned int > component_to_base_index (const unsigned int component) const
 
std::pair< unsigned int, unsigned int > block_to_base_index (const unsigned int block) const
 
std::pair< unsigned int, types::global_dof_indexsystem_to_block_index (const unsigned int component) const
 
unsigned int component_to_block_index (const unsigned int component) const
 
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
 
const std::vector< Point< dim > > & get_unit_support_points () const
 
bool has_support_points () const
 
virtual Point< dim > unit_support_point (const unsigned int index) const
 
const std::vector< Point< dim-1 > > & get_unit_face_support_points () const
 
bool has_face_support_points () const
 
virtual Point< dim-1 > unit_face_support_point (const unsigned int index) const
 
const std::vector< Point< dim > > & get_generalized_support_points () const
 
bool has_generalized_support_points () const
 
const std::vector< Point< dim-1 > > & get_generalized_face_support_points () const
 
bool has_generalized_face_support_points () const
 
GeometryPrimitive get_associated_geometry_primitive (const unsigned int cell_dof_index) const
 
- Public Member Functions inherited from Subscriptor
 Subscriptor ()
 
 Subscriptor (const Subscriptor &)
 
 Subscriptor (Subscriptor &&) noexcept
 
virtual ~Subscriptor ()
 
Subscriptoroperator= (const Subscriptor &)
 
Subscriptoroperator= (Subscriptor &&) noexcept
 
void subscribe (const char *identifier=nullptr) const
 
void unsubscribe (const char *identifier=nullptr) const
 
unsigned int n_subscriptions () const
 
void list_subscribers () const
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 
- Public Member Functions inherited from FiniteElementData< dim >
 FiniteElementData (const std::vector< unsigned int > &dofs_per_object, const unsigned int n_components, const unsigned int degree, const Conformity conformity=unknown, const BlockIndices &block_indices=BlockIndices())
 
unsigned int n_dofs_per_vertex () const
 
unsigned int n_dofs_per_line () const
 
unsigned int n_dofs_per_quad () const
 
unsigned int n_dofs_per_hex () const
 
unsigned int n_dofs_per_face () const
 
unsigned int n_dofs_per_cell () const
 
template<int structdim>
unsigned int n_dofs_per_object () const
 
unsigned int n_components () const
 
unsigned int n_blocks () const
 
const BlockIndicesblock_indices () const
 
unsigned int tensor_degree () const
 
bool conforms (const Conformity) const
 
bool operator== (const FiniteElementData &) const
 

Static Private Member Functions

static std::vector< bool > get_riaf_vector (const unsigned int degree)
 
static std::vector< unsigned int > get_dpo_vector (const unsigned int degree)
 

Private Attributes

const unsigned int n_bubbles
 

Additional Inherited Members

- Public Types inherited from FiniteElementData< dim >
enum  Conformity {
  unknown = 0x00, L2 = 0x01, Hcurl = 0x02, Hdiv = 0x04,
  H1 = Hcurl | Hdiv, H2 = 0x0e
}
 
- Static Public Member Functions inherited from FE_Q_Base< TensorProductPolynomialsBubbles< dim >, dim, spacedim >
static ::ExceptionBaseExcFEQCannotHaveDegree0 ()
 
- Static Public Member Functions inherited from FiniteElement< dim, spacedim >
static ::ExceptionBaseExcShapeFunctionNotPrimitive (int arg1)
 
static ::ExceptionBaseExcFENotPrimitive ()
 
static ::ExceptionBaseExcUnitShapeValuesDoNotExist ()
 
static ::ExceptionBaseExcFEHasNoSupportPoints ()
 
static ::ExceptionBaseExcEmbeddingVoid ()
 
static ::ExceptionBaseExcProjectionVoid ()
 
static ::ExceptionBaseExcWrongInterfaceMatrixSize (int arg1, int arg2)
 
static ::ExceptionBaseExcInterpolationNotImplemented ()
 
- Static Public Member Functions inherited from Subscriptor
static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 
- Public Attributes inherited from FiniteElementData< dim >
const unsigned int dofs_per_vertex
 
const unsigned int dofs_per_line
 
const unsigned int dofs_per_quad
 
const unsigned int dofs_per_hex
 
const unsigned int first_line_index
 
const unsigned int first_quad_index
 
const unsigned int first_hex_index
 
const unsigned int first_face_line_index
 
const unsigned int first_face_quad_index
 
const unsigned int dofs_per_face
 
const unsigned int dofs_per_cell
 
const unsigned int components
 
const unsigned int degree
 
const Conformity conforming_space
 
const BlockIndices block_indices_data
 
- Static Public Attributes inherited from FiniteElement< dim, spacedim >
static const unsigned int space_dimension = spacedim
 
- Static Public Attributes inherited from FiniteElementData< dim >
static const unsigned int dimension = dim
 
- Protected Member Functions inherited from FE_Q_Base< TensorProductPolynomialsBubbles< dim >, dim, spacedim >
void initialize (const std::vector< Point< 1 > > &support_points_1d)
 
void initialize_constraints (const std::vector< Point< 1 > > &points)
 
void initialize_unit_support_points (const std::vector< Point< 1 > > &points)
 
void initialize_unit_face_support_points (const std::vector< Point< 1 > > &points)
 
void initialize_quad_dof_index_permutation ()
 
- Protected Member Functions inherited from FE_Poly< TensorProductPolynomialsBubbles< dim >, dim, spacedim >
void correct_third_derivatives (internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data, const internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const unsigned int n_q_points, const unsigned int dof) const
 
- Protected Member Functions inherited from FiniteElement< dim, spacedim >
void reinit_restriction_and_prolongation_matrices (const bool isotropic_restriction_only=false, const bool isotropic_prolongation_only=false)
 
TableIndices< 2 > interface_constraints_size () const
 
virtual std::unique_ptr< InternalDataBaseget_face_data (const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim-1 > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const
 
virtual std::unique_ptr< InternalDataBaseget_subface_data (const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim-1 > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const
 
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 InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const =0
 
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 InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const =0
 
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 InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const =0
 
- Static Protected Member Functions inherited from FE_Q_Base< TensorProductPolynomialsBubbles< dim >, dim, spacedim >
static std::vector< unsigned int > get_dpo_vector (const unsigned int degree)
 
- Static Protected Member Functions inherited from FiniteElement< dim, spacedim >
static std::vector< unsigned int > compute_n_nonzero_components (const std::vector< ComponentMask > &nonzero_components)
 
- Protected Attributes inherited from FE_Poly< TensorProductPolynomialsBubbles< dim >, dim, spacedim >
TensorProductPolynomialsBubbles< dim > poly_space
 
- Protected Attributes inherited from FiniteElement< dim, spacedim >
std::vector< std::vector< FullMatrix< double > > > restriction
 
std::vector< std::vector< FullMatrix< double > > > prolongation
 
FullMatrix< double > interface_constraints
 
std::vector< Point< dim > > unit_support_points
 
std::vector< Point< dim-1 > > unit_face_support_points
 
std::vector< Point< dim > > generalized_support_points
 
std::vector< Point< dim-1 > > generalized_face_support_points
 
Table< 2, int > adjust_quad_dof_index_for_face_orientation_table
 
std::vector< int > adjust_line_dof_index_for_line_orientation_table
 
std::vector< std::pair< unsigned int, unsigned int > > system_to_component_table
 
std::vector< std::pair< unsigned int, unsigned int > > face_system_to_component_table
 
std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > system_to_base_table
 
std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > face_system_to_base_table
 
BlockIndices base_to_block_indices
 
std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > component_to_base_table
 
const std::vector< bool > restriction_is_additive_flags
 
const std::vector< ComponentMasknonzero_components
 
const std::vector< unsigned int > n_nonzero_components_table
 
const bool cached_primitivity
 

Detailed Description

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

Implementation of a scalar Lagrange finite element Q_p^+ that yields the finite element space of continuous, piecewise polynomials of degree p in each coordinate direction plus some bubble enrichment space spanned by \((2x_j-1)^{p-1}\prod_{i=0}^{dim-1}(x_i(1-x_i))\). Therefore the highest polynomial degree is \(p+1\). This class is realized using tensor product polynomials based on equidistant or given support points.

The standard constructor of this class takes the degree p of this finite element. Alternatively, it can take a quadrature formula points defining the support points of the Lagrange interpolation in one coordinate direction.

For more information about the spacedim template parameter check the documentation of FiniteElement or the one of Triangulation.

Due to the fact that the enrichments are small almost everywhere for large p, the condition number for the mass and stiffness matrix fastly increaseses with increasing p. Below you see a comparison with FE_Q(QGaussLobatto(p+1)) for dim=1.

fe_q_bubbles_conditioning.png

Therefore, this element should be used with care for \(p>3\).

Implementation

The constructor creates a TensorProductPolynomials object that includes the tensor product of LagrangeEquidistant polynomials of degree p plus the bubble enrichments. This TensorProductPolynomialsBubbles object provides all values and derivatives of the shape functions. In case a quadrature rule is given, the constructor creates a TensorProductPolynomialsBubbles object that includes the tensor product of Lagrange polynomials with the support points from points and the bubble enrichments as defined above.

Furthermore the constructor fills the interface_constrains, the prolongation (embedding) and the restriction matrices.

Numbering of the degrees of freedom (DoFs)

The original ordering of the shape functions represented by the TensorProductPolynomialsBubbles is a tensor product numbering. However, the shape functions on a cell are renumbered beginning with the shape functions whose support points are at the vertices, then on the line, on the quads, and finally (for 3d) on the hexes. Finally, there are support points for the bubble enrichments in the middle of the cell.

Definition at line 82 of file fe_q_bubbles.h.

Constructor & Destructor Documentation

◆ FE_Q_Bubbles() [1/2]

template<int dim, int spacedim>
FE_Q_Bubbles< dim, spacedim >::FE_Q_Bubbles ( const unsigned int  p)

Constructor for tensor product polynomials of degree p plus bubble enrichments

Definition at line 176 of file fe_q_bubbles.cc.

◆ FE_Q_Bubbles() [2/2]

template<int dim, int spacedim>
FE_Q_Bubbles< dim, spacedim >::FE_Q_Bubbles ( const Quadrature< 1 > &  points)

Constructor for tensor product polynomials with support points points plus bubble enrichments based on a one-dimensional quadrature formula. The degree of the finite element is points.size(). Note that the first point has to be 0 and the last one 1.

Definition at line 214 of file fe_q_bubbles.cc.

Member Function Documentation

◆ get_name()

template<int dim, int spacedim>
std::string FE_Q_Bubbles< dim, spacedim >::get_name ( ) const
virtual

Return a string that uniquely identifies a finite element. This class returns FE_Q_Bubbles<dim>(degree), with dim and degree replaced by appropriate values.

Implements FiniteElement< dim, spacedim >.

Definition at line 252 of file fe_q_bubbles.cc.

◆ convert_generalized_support_point_values_to_dof_values()

template<int dim, int spacedim>
void FE_Q_Bubbles< dim, spacedim >::convert_generalized_support_point_values_to_dof_values ( const std::vector< Vector< double > > &  support_point_values,
std::vector< double > &  nodal_values 
) const
virtual

Given the values of a function \(f(\mathbf x)\) at the (generalized) support points of the reference cell, this function then computes what the nodal values of the element are, i.e., \(\Psi_i[f]\), where \(\Psi_i\) are the node functionals of the element (see also Node values or node functionals). The values \(\Psi_i[f]\) are then the expansion coefficients for the shape functions of the finite element function that interpolates the given function \(f(x)\), i.e., \( f_h(\mathbf x) = \sum_i \Psi_i[f] \varphi_i(\mathbf x) \) is the finite element interpolant of \(f\) with the current element. The operation described here is used, for example, in the FETools::compute_node_matrix() function.

In more detail, let us assume that the generalized support points (see this glossary entry ) of the current element are \(\hat{\mathbf x}_i\) and that the node functionals associated with the current element are \(\Psi_i[\cdot]\). Then, the fact that the element is based on generalized support points, implies that if we apply \(\Psi_i\) to a (possibly vector-valued) finite element function \(\varphi\), the result must have the form \(\Psi_i[\varphi] = f_i(\varphi(\hat{\mathbf x}_i))\) – in other words, the value of the node functional \(\Psi_i\) applied to \(\varphi\) only depends on the values of \(\varphi\) at \(\hat{\mathbf x}_i\) and not on values anywhere else, or integrals of \(\varphi\), or any other kind of information.

The exact form of \(f_i\) depends on the element. For example, for scalar Lagrange elements, we have that in fact \(\Psi_i[\varphi] = \varphi(\hat{\mathbf x}_i)\). If you combine multiple scalar Lagrange elements via an FESystem object, then \(\Psi_i[\varphi] = \varphi(\hat{\mathbf x}_i)_{c(i)}\) where \(c(i)\) is the result of the FiniteElement::system_to_component_index() function's return value's first component. In these two cases, \(f_i\) is therefore simply the identity (in the scalar case) or a function that selects a particular vector component of its argument. On the other hand, for Raviart-Thomas elements, one would have that \(f_i(\mathbf y) = \mathbf y \cdot \mathbf n_i\) where \(\mathbf n_i\) is the normal vector of the face at which the shape function is defined.

Given all of this, what this function does is the following: If you input a list of values of a function \(\varphi\) at all generalized support points (where each value is in fact a vector of values with as many components as the element has), then this function returns a vector of values obtained by applying the node functionals to these values. In other words, if you pass in \(\{\varphi(\hat{\mathbf x}_i)\}_{i=0}^{N-1}\) then you will get out a vector \(\{\Psi[\varphi]\}_{i=0}^{N-1}\) where \(N\) equals dofs_per_cell.

Parameters
[in]support_point_valuesAn array of size dofs_per_cell (which equals the number of points the get_generalized_support_points() function will return) where each element is a vector with as many entries as the element has vector components. This array should contain the values of a function at the generalized support points of the current element.
[out]nodal_valuesAn array of size dofs_per_cell that contains the node functionals of the element applied to the given function.
Note
It is safe to call this function for (transformed) values on the real cell only for elements with trivial MappingType. For all other elements (for example for H(curl), or H(div) conforming elements) vector values have to be transformed to the reference cell first.
Given what the function is supposed to do, the function clearly can only work for elements that actually implement (generalized) support points. Elements that do not have generalized support points – e.g., elements whose nodal functionals evaluate integrals or moments of functions (such as FE_Q_Hierarchical) – can in general not make sense of the operation that is required for this function. They consequently may not implement it.

Reimplemented from FiniteElement< dim, spacedim >.

Definition at line 339 of file fe_q_bubbles.cc.

◆ get_interpolation_matrix()

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

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 only available if the source element is also a FE_Q_Bubbles element. Otherwise, an exception of type FiniteElement<dim,spacedim>::ExcInterpolationNotImplemented is thrown.

Reimplemented from FE_Q_Base< TensorProductPolynomialsBubbles< dim >, dim, spacedim >.

Definition at line 367 of file fe_q_bubbles.cc.

◆ get_prolongation_matrix()

template<int dim, int spacedim>
const FullMatrix< double > & FE_Q_Bubbles< dim, spacedim >::get_prolongation_matrix ( const unsigned int  child,
const RefinementCase< dim > &  refinement_case 
) const
virtual

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 projection matrices are not implemented in the derived finite element class, this function aborts with ExcEmbeddingVoid. You can check whether this is the case by calling the prolongation_is_implemented() or the isotropic_prolongation_is_implemented() function.

Reimplemented from FE_Q_Base< TensorProductPolynomialsBubbles< dim >, dim, spacedim >.

Definition at line 444 of file fe_q_bubbles.cc.

◆ get_restriction_matrix()

template<int dim, int spacedim>
const FullMatrix< double > & FE_Q_Bubbles< dim, spacedim >::get_restriction_matrix ( const unsigned int  child,
const RefinementCase< dim > &  refinement_case 
) const
virtual

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 ExcProjectionVoid. You can check whether this is the case by calling the restriction_is_implemented() or the isotropic_restriction_is_implemented() function.

Reimplemented from FE_Q_Base< TensorProductPolynomialsBubbles< dim >, dim, spacedim >.

Definition at line 465 of file fe_q_bubbles.cc.

◆ has_support_on_face()

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

Check for non-zero values on a face.

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

Implementation of the interface in FiniteElement

Reimplemented from FE_Q_Base< TensorProductPolynomialsBubbles< dim >, dim, spacedim >.

Definition at line 429 of file fe_q_bubbles.cc.

◆ clone()

template<int dim, int spacedim>
std::unique_ptr< FiniteElement< dim, spacedim > > FE_Q_Bubbles< dim, spacedim >::clone ( ) const
virtual

A sort of virtual copy constructor, this function returns a copy of the finite element object. Derived classes need to override the function here in this base class and return an object of the same type as the derived class.

Some places in the library, for example the constructors of FESystem as well as the hp::FECollection class, need to make copies of finite elements without knowing their exact type. They do so through this function.

Implements FiniteElement< dim, spacedim >.

Definition at line 329 of file fe_q_bubbles.cc.

◆ get_riaf_vector()

template<int dim, int spacedim>
std::vector< bool > FE_Q_Bubbles< dim, spacedim >::get_riaf_vector ( const unsigned int  degree)
staticprivate

Return the restriction_is_additive flags. Only the last components for the bubble enrichments are true.

Definition at line 403 of file fe_q_bubbles.cc.

◆ get_dpo_vector()

template<int dim, int spacedim>
std::vector< unsigned int > FE_Q_Bubbles< dim, spacedim >::get_dpo_vector ( const unsigned int  degree)
staticprivate

Only for internal use. Its full name is get_dofs_per_object_vector function and it creates the dofs_per_object vector that is needed within the constructor to be passed to the constructor of FiniteElementData.

Definition at line 414 of file fe_q_bubbles.cc.

Member Data Documentation

◆ n_bubbles

template<int dim, int spacedim = dim>
const unsigned int FE_Q_Bubbles< dim, spacedim >::n_bubbles
private

Number of additional bubble functions

Definition at line 168 of file fe_q_bubbles.h.


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