Reference documentation for deal.II version 9.1.1
|
#include <deal.II/fe/fe_q_iso_q1.h>
Public Member Functions | |
FE_Q_iso_Q1 (const unsigned int n_subdivisions) | |
virtual std::string | get_name () const override |
virtual std::unique_ptr< FiniteElement< dim, spacedim > > | clone () 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 |
Functions to support hp | |
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_Q_Base< TensorProductPolynomials< dim, Polynomials::PiecewisePolynomial< double > >, dim, spacedim > | |
FE_Q_Base (const TensorProductPolynomials< dim, Polynomials::PiecewisePolynomial< double > > &poly_space, const FiniteElementData< dim > &fe_data, const std::vector< bool > &restriction_is_additive_flags) | |
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 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 std::pair< Table< 2, bool >, std::vector< unsigned int > > | get_constant_modes () const override |
virtual bool | hp_constraints_are_implemented () 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 |
Public Member Functions inherited from FE_Poly< TensorProductPolynomials< dim, Polynomials::PiecewisePolynomial< double > >, dim, spacedim > | |
FE_Poly (const TensorProductPolynomials< dim, Polynomials::PiecewisePolynomial< double > > &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 |
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 |
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 |
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 |
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< TensorProductPolynomials< dim, Polynomials::PiecewisePolynomial< double > >, dim, spacedim > | |
static ::ExceptionBase & | ExcFEQCannotHaveDegree0 () |
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 |
Protected Member Functions inherited from FE_Q_Base< TensorProductPolynomials< dim, Polynomials::PiecewisePolynomial< double > >, 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< TensorProductPolynomials< dim, Polynomials::PiecewisePolynomial< double > >, 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 |
Static Protected Member Functions inherited from FE_Q_Base< TensorProductPolynomials< dim, Polynomials::PiecewisePolynomial< double > >, 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< TensorProductPolynomials< dim, Polynomials::PiecewisePolynomial< double > >, dim, spacedim > | |
TensorProductPolynomials< dim, Polynomials::PiecewisePolynomial< double > > | 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 a scalar Lagrange finite element Qp-iso-Q1
that defines the finite element space of continuous, piecewise linear elements with p
subdivisions in each coordinate direction. It yields an element with the same number of degrees of freedom as the Qp
elements but using linear interpolation instead of higher order one. This type of element is also called macro element in the literature as it really consists of several smaller elements, namely pdim
such sub-cells.
The numbering of degrees of freedom is done in exactly the same way as in FE_Q of degree p
. See there for a detailed description on how degrees of freedom are numbered within one element.
This element represents a Q-linear finite element space on a reduced mesh size h/p. Its effect is equivalent to using FE_Q of degree one on a finer mesh by a factor p
if an equivalent quadrature is used. However, this element reduces the flexibility in the choice of (adaptive) mesh size by exactly this factor p
, which typically reduces efficiency. On the other hand, comparing this element with p
subdivisions to the FE_Q element of degree p
on the same mesh shows that the convergence is typically much worse for smooth problems. In particular, Qp
elements achieve interpolation orders of hp+1 in the L2 norm, whereas these elements reach only (h/p)2. For these two reasons, this element is usually not very useful as a standalone. In addition, any evaluation of face terms on the boundaries within the elements becomes impossible with this element because deal.II does not have the equivalent of FEFaceValues for lower-dimensional integrals in the interior of cells.
Nonetheless, there are a few use cases where this element actually is useful:
Systems of PDEs where certain variables demand for higher resolutions than the others and the additional degrees of freedom should be spent on increasing the resolution of linears instead of higher order polynomials, and you do not want to use two different meshes for the different components. This can be the case when irregularities (shocks) appear in the solution and stabilization techniques are used that work for linears but not higher order elements.
Stokes/Navier Stokes systems such as the one discussed in step-22 could be solved with Q2-iso-Q1 elements for velocities instead of Q2 elements. Combined with Q1 pressures they give a stable mixed element pair. However, they perform worse than the standard (Taylor-Hood \(Q_2\times Q_1\)) approach in most situations.
p
with a preconditioner based on Qp-iso-Q1
elements: Some preconditioners like algebraic multigrid perform much better with linear elements than with higher order elements because they often implicitly assume a sparse connectivity between entries. Then, creating a preconditioner matrix based on these elements yields the same number of degrees of freedom (and a spectrally equivalent linear system), which can be combined with a (high order) system matrix in an iterative solver like CG. Due to the nature of these elements as a concatenation of linears, care must be taken when selecting quadrature formulas for this element. The standard choice for an element of p
subelements is a formula QIterated<dim>(QGauss<1>(2), p)
, which corresponds to the formula that would be used for integrating functions on a finer mesh. This is in contrast with FE_Q(p) where QGauss<dim>(p+1) is the default choice. In particular, care must be taken to not use a quadrature formula that evaluates the basis functions (and their derivatives) on sub-element boundaries as the gradients of piecewiese functions on internal boundaries are set to zero. No checks are performed internally to ensure that this is not the case - it is the user's responsibility to avoid these situations.
Also note that the usual deal.II routines for setting up sparsity patterns and assembling matrices do not make use of the increased sparsity in this element compared to FE_Q. This is because DoFTools::make_sparsity_pattern assumes coupling between all degrees of freedom within the element, whereas FE_Q_iso_Q1 with more than one subdivision does have less coupling.
Definition at line 113 of file fe_q_iso_q1.h.
FE_Q_iso_Q1< dim, spacedim >::FE_Q_iso_Q1 | ( | const unsigned int | n_subdivisions | ) |
Construct a FE_Q_iso_Q1 element with a given number of subdivisions. The number of subdivision is similar to the degree in FE_Q in the sense that both elements produce the same number of degrees of freedom.
Definition at line 34 of file fe_q_iso_q1.cc.
|
overridevirtual |
Return a string that uniquely identifies a finite element. This class returns FE_Q_iso_q1<dim>(equivalent_degree)
, with dim
and equivalent_degree
replaced by appropriate values.
Implements FiniteElement< dim, spacedim >.
Definition at line 63 of file fe_q_iso_q1.cc.
|
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 >.
Definition at line 101 of file fe_q_iso_q1.cc.
|
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 >.
Definition at line 80 of file fe_q_iso_q1.cc.
|
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 >.
Definition at line 110 of file fe_q_iso_q1.cc.