Reference documentation for deal.II version 9.1.1
|
#include <deal.II/fe/fe_dgq.h>
Public Member Functions | |
FE_DGQ (const unsigned int p) | |
virtual std::string | get_name () const override |
virtual void | get_interpolation_matrix (const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) 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 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 bool | has_support_on_face (const unsigned int shape_index, const unsigned int face_index) const override |
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > | get_constant_modes () const override |
virtual void | convert_generalized_support_point_values_to_dof_values (const std::vector< Vector< double >> &support_point_values, std::vector< double > &nodal_values) const override |
virtual std::size_t | memory_consumption () const override |
virtual std::unique_ptr< FiniteElement< dim, spacedim > > | clone () const override |
Functions to support hp | |
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 bool | hp_constraints_are_implemented () const override |
virtual FiniteElementDomination::Domination | compare_for_domination (const FiniteElement< dim, spacedim > &fe_other, const unsigned int codim=0) const override final |
Public Member Functions inherited from FE_Poly< TensorProductPolynomials< dim >, dim, spacedim > | |
FE_Poly (const TensorProductPolynomials< 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 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 |
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 () override=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 |
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 |
virtual FiniteElementDomination::Domination | compare_for_face_domination (const FiniteElement< dim, spacedim > &fe_other) const final |
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 |
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 |
unsigned int | adjust_line_dof_index_for_line_orientation (const unsigned int index, const bool line_orientation) const |
const ComponentMask & | get_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_index > | system_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 () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (Subscriptor &&) noexcept |
void | subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
void | unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
unsigned int | n_subscriptions () const |
template<typename StreamType > | |
void | list_subscribers (StreamType &stream) const |
void | list_subscribers () const |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
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 BlockIndices & | block_indices () const |
unsigned int | tensor_degree () const |
bool | conforms (const Conformity) const |
bool | operator== (const FiniteElementData &) const |
Protected Member Functions | |
FE_DGQ (const std::vector< Polynomials::Polynomial< double >> &polynomials) | |
Protected Member Functions inherited from FE_Poly< TensorProductPolynomials< 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< 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 |
virtual std::unique_ptr< 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 |
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 |
Private Member Functions | |
void | rotate_indices (std::vector< unsigned int > &indices, const char direction) const |
Static Private Member Functions | |
static std::vector< unsigned int > | get_dpo_vector (const unsigned int degree) |
Friends | |
template<int dim1, int spacedim1> | |
class | FE_DGQ |
template<int dim1, int spacedim1> | |
class | MappingQ |
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 FiniteElement< dim, spacedim > | |
static ::ExceptionBase & | ExcShapeFunctionNotPrimitive (int arg1) |
static ::ExceptionBase & | ExcFENotPrimitive () |
static ::ExceptionBase & | ExcUnitShapeValuesDoNotExist () |
static ::ExceptionBase & | ExcFEHasNoSupportPoints () |
static ::ExceptionBase & | ExcEmbeddingVoid () |
static ::ExceptionBase & | ExcProjectionVoid () |
static ::ExceptionBase & | ExcWrongInterfaceMatrixSize (int arg1, int arg2) |
static ::ExceptionBase & | ExcInterpolationNotImplemented () |
Static Public Member Functions inherited from Subscriptor | |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (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 |
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< TensorProductPolynomials< dim >, dim, spacedim > | |
TensorProductPolynomials< 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< ComponentMask > | nonzero_components |
const std::vector< unsigned int > | n_nonzero_components_table |
const bool | cached_primitivity |
Implementation of scalar, discontinuous tensor product elements based on equidistant support points.
This is a discontinuous finite element based on tensor products of Lagrangian polynomials. The shape functions are Lagrangian interpolants of an equidistant grid of points on the unit cell. The points are numbered in lexicographical order, with x running fastest, then y, then z (if these coordinates are present for a given space dimension at all). For example, these are the node orderings for FE_DGQ(1)
in 3d:
* 6-------7 6-------7 * /| | / /| * / | | / / | * / | | / / | * 4 | | 4-------5 | * | 2-------3 | | 3 * | / / | | / * | / / | | / * |/ / | |/ * 0-------1 0-------1 *
and FE_DGQ(2)
:
* 24--25--26 24--25--26 * /| | / /| * 21 | | 21 22 23 | * / 15 16 17 / / 17 * 18 | | 18--19--20 | * |12 6---7---8 | |14 8 * 9 / / 9 10 11 / * | 3 4 5 | | 5 * |/ / | |/ * 0---1---2 0---1---2 *
with node 13 being placed in the interior of the hex.
Note, however, that these are just the Lagrange interpolation points of the shape functions. Even though they may physically be on the boundary of the cell, they are logically in the interior since there are no continuity requirements for these shape functions across cell boundaries. While discontinuous, when restricted to a single cell the shape functions of this element are exactly the same as those of the FE_Q element where they are shown visually.
When constructing an FE_DGQ element at polynomial degrees one or two, equidistant support points at 0 and 1 (linear case) or 0, 0.5, and 1 (quadratic case) are used. The unit support or nodal points xi are those points where the jth Lagrange polynomial satisfies the \(\delta_{ij}\) property, i.e., where one polynomial is one and all the others are zero. For higher polynomial degrees, the support points are non-equidistant by default, and chosen to be the support points of the (degree+1)
-order Gauss-Lobatto quadrature rule. This point distribution yields well-conditioned Lagrange interpolation at arbitrary polynomial degrees. By contrast, polynomials based on equidistant points get increasingly ill-conditioned as the polynomial degree increases. In interpolation, this effect is known as the Runge phenomenon. For Galerkin methods, the Runge phenomenon is typically not visible in the solution quality but rather in the condition number of the associated system matrices. For example, the elemental mass matrix of equidistant points at degree 10 has condition number 2.6e6, whereas the condition number for Gauss-Lobatto points is around 400.
The Gauss-Lobatto points in 1D include the end points 0 and +1 of the unit interval. The interior points are shifted towards the end points, which gives a denser point distribution close to the element boundary.
|
protected |
Constructor for tensor product polynomials based on an arbitrary vector of polynomials. This constructor is used in derived classes to construct e.g. elements with arbitrary nodes or elements based on Legendre polynomials.
The degree of these polynomials is polynomials.size()-1
.
|
overridevirtual |
Return a string that uniquely identifies a finite element. This class returns FE_DGQ<dim>(degree)
, with dim
and degree
replaced by appropriate values.
Implements FiniteElement< dim, spacedim >.
Reimplemented in FE_DGQHermite< dim, spacedim >, FE_DGQLegendre< dim, spacedim >, and FE_DGQArbitraryNodes< dim, spacedim >.
|
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 only available if the source element is also a FE_DGQ
element. Otherwise, an exception of type FiniteElement<dim>::ExcInterpolationNotImplemented is thrown.
Reimplemented from FiniteElement< dim, spacedim >.
|
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
.
Derived elements 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>::ExcInterpolationNotImplemented.
Reimplemented from FiniteElement< dim, spacedim >.
|
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
.
Derived elements 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>::ExcInterpolationNotImplemented.
Reimplemented from FiniteElement< dim, spacedim >.
|
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.
Reimplemented from FiniteElement< dim, spacedim >.
|
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.
Reimplemented from FiniteElement< dim, spacedim >.
|
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.
This being a discontinuous element, the set of such constraints is of course empty.
Reimplemented from FiniteElement< dim, spacedim >.
|
overridevirtual |
Same as hp_vertex_dof_indices(), except that the function treats degrees of freedom on lines.
This being a discontinuous element, the set of such constraints is of course empty.
Reimplemented from FiniteElement< dim, spacedim >.
|
overridevirtual |
Same as hp_vertex_dof_indices(), except that the function treats degrees of freedom on quads.
This being a discontinuous element, the set of such constraints is of course empty.
Reimplemented from FiniteElement< dim, spacedim >.
|
overridevirtual |
Return whether this element implements its hanging node constraints in the new way, which has to be used to make elements "hp compatible".
For the FE_DGQ class the result is always true (independent of the degree of the element), as it has no hanging nodes (being a discontinuous element).
Reimplemented from FiniteElement< dim, spacedim >.
|
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.
Reimplemented from FiniteElement< dim, spacedim >.
|
overridevirtual |
This function returns true
, if the shape function shape_index
has non-zero function values somewhere on the face face_index
.
Reimplemented from FiniteElement< dim, spacedim >.
|
overridevirtual |
Return a list of constant modes of the element. For this element, it simply returns one row with all entries set to true.
Reimplemented from FiniteElement< dim, spacedim >.
Reimplemented in FE_DGQLegendre< dim, spacedim >.
|
overridevirtual |
Implementation of the corresponding function in the FiniteElement class. Since the current element is interpolatory, the nodal values are exactly the support point values. Furthermore, since the current element is scalar, the support point values need to be vectors of length 1.
Reimplemented from FiniteElement< dim, spacedim >.
Reimplemented in FE_DGQArbitraryNodes< dim, spacedim >.
|
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.
Reimplemented from FiniteElement< dim, spacedim >.
|
overridevirtual |
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 >.
Reimplemented in FE_DGQHermite< dim, spacedim >, FE_DGQLegendre< dim, spacedim >, and FE_DGQArbitraryNodes< dim, spacedim >.
|
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
.
|
private |
Compute renumbering for rotation of degrees of freedom.
This function rotates a tensor product numbering of degrees of freedom by 90 degrees. It is used to compute the transfer matrices of the children by using only the matrix for the first child.
The direction parameter determines the type of rotation. It is one character of xXyYzZ
. The character determines the axis of rotation, case determines the direction. Lower case is counter-clockwise seen in direction of the axis.
Since rotation around the y-axis is not used, it is not implemented either.
|
friend |