Reference documentation for deal.II version 9.1.1
|
#include <deal.II/fe/fe_enriched.h>
Classes | |
class | InternalData |
Public Member Functions | |
FE_Enriched (const FiniteElement< dim, spacedim > &fe_base, const FiniteElement< dim, spacedim > &fe_enriched, const Function< spacedim > *enrichment_function) | |
FE_Enriched (const FiniteElement< dim, spacedim > &fe_base) | |
FE_Enriched (const FiniteElement< dim, spacedim > *fe_base, const std::vector< const FiniteElement< dim, spacedim > *> &fe_enriched, const std::vector< std::vector< std::function< const Function< spacedim > *(const typename Triangulation< dim, spacedim >::cell_iterator &)>>> &functions) | |
virtual std::unique_ptr< FiniteElement< dim, spacedim > > | clone () const override |
virtual UpdateFlags | requires_update_flags (const UpdateFlags update_flags) const override |
virtual std::string | get_name () const override |
virtual const FiniteElement< dim, spacedim > & | base_element (const unsigned int index) const override |
virtual double | shape_value (const unsigned int i, const Point< dim > &p) const override |
const std::vector< std::vector< std::function< const Function< spacedim > *(const typename Triangulation< dim, spacedim >::cell_iterator &)> > > | get_enrichments () const |
const FESystem< dim, spacedim > & | get_fe_system () const |
Transfer matrices | |
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 |
Functions to support hp | |
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 |
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 |
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 |
virtual bool | has_support_on_face (const unsigned int shape_index, const unsigned int face_index) 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 void | get_interpolation_matrix (const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) 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 |
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 |
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > | get_constant_modes () 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 |
virtual void | convert_generalized_support_point_values_to_dof_values (const std::vector< Vector< double >> &support_point_values, std::vector< double > &nodal_values) 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 | |
template<int dim_1> | |
std::unique_ptr< typename FiniteElement< dim, spacedim >::InternalDataBase > | setup_data (std::unique_ptr< typename FESystem< dim, spacedim >::InternalData > fes_data, const UpdateFlags flags, const Quadrature< dim_1 > &quadrature) const |
virtual std::unique_ptr< typename FiniteElement< dim, spacedim >::InternalDataBase > | get_data (const UpdateFlags 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 |
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 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 |
Protected Attributes | |
std::vector< std::vector< std::vector< unsigned int > > > | base_no_mult_local_enriched_dofs |
const std::vector< std::vector< std::function< const Function< spacedim > *(const typename Triangulation< dim, spacedim >::cell_iterator &)> > > | enrichments |
const bool | is_enriched |
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 |
Private Member Functions | |
FE_Enriched (const std::vector< const FiniteElement< dim, spacedim > *> &fes, const std::vector< unsigned int > &multiplicities, const std::vector< std::vector< std::function< const Function< spacedim > *(const typename Triangulation< dim, spacedim >::cell_iterator &)>>> &functions) | |
void | initialize (const std::vector< const FiniteElement< dim, spacedim > *> &fes, const std::vector< unsigned int > &multiplicities) |
template<int dim_1> | |
void | multiply_by_enrichment (const Quadrature< dim_1 > &quadrature, const InternalData &fe_data, const internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const typename Triangulation< dim, spacedim >::cell_iterator &cell, internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const |
Private Attributes | |
const std::unique_ptr< const FESystem< dim, spacedim > > | fe_system |
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) |
Implementation of a partition of unity finite element method (PUM) by Babuska and Melenk which enriches a standard finite element with an enrichment function multiplied with another (usually linear) finite element:
\[ U(\mathbf x) = \sum_i N_i(\mathbf x) U_i + \sum_j N_j(\mathbf x) \sum_k F_k(\mathbf x) U_{jk} \]
where \( N_i(\mathbf x) \) and \( N_j(\mathbf x) \) are the underlying finite elements (including the mapping from the isoparametric element to the real element); \( F_k(\mathbf x) \) are the scalar enrichment functions in real space (e.g. \( 1/r \), \( \exp(-r) \), etc); \( U_i \) and \( U_{jk} \) are the standard and enriched DoFs. This allows to include in the finite element space a priori knowledge about the partial differential equation being solved which in turn improves the local approximation properties of the spaces. This can be useful for highly oscillatory solutions, problems with domain corners or on unbounded domains or sudden changes of boundary conditions. PUM method uses finite element spaces which satisfy the partition of unity property (e.g. FE_Q). Among other properties this makes the resulting space to reproduce enrichment functions exactly.
The simplest constructor of this class takes two finite element objects and an enrichment function to be used. For example
In this case, standard DoFs are distributed by FE_Q<dim>(2)
, whereas enriched DoFs are coming from a single finite element FE_Q<dim>(1)
used with a single enrichment function function
. In this case, the total number of DoFs on the enriched element is the sum of DoFs from FE_Q<dim>(2)
and FE_Q<dim>(1)
.
As an example of an enrichment function, consider \( \exp(-x) \), which leads to the following shape functions on the unit element:
1d element, base and enriched shape functions. | enriched shape function corresponding to the central vertex. |
Note that evaluation of gradients (hessians) of the enriched shape functions or the finite element field requires evaluation of gradients (gradients and hessians) of the enrichment functions:
\begin{align*} U(\mathbf x) &= \sum_i N_i(\mathbf x) U_i + \sum_{j,k} N_j(\mathbf x) F_k(\mathbf x) U_{jk} \\ \mathbf \nabla U(\mathbf x) &= \sum_i \mathbf \nabla N_i(\mathbf x) U_i + \sum_{j,k} \left[\mathbf \nabla N_j(\mathbf x) F_k(\mathbf x) + N_j(\mathbf x) \mathbf \nabla F_k(\mathbf x) \right] U_{jk} \\ \mathbf \nabla \mathbf \nabla U(\mathbf x) &= \sum_i \mathbf \nabla \mathbf \nabla N_i(\mathbf x) U_i + \sum_{j,k} \left[\mathbf \nabla \mathbf \nabla N_j(\mathbf x) F_k(\mathbf x) + \mathbf \nabla F_k(\mathbf x) \mathbf \nabla N_j(\mathbf x) + \mathbf \nabla N_j(\mathbf x) \mathbf \nabla F_k(\mathbf x) + N_j(\mathbf x) \mathbf \nabla \mathbf \nabla F_k(\mathbf x) \right] U_{jk} \end{align*}
In most applications it is beneficial to introduce enrichments only in some part of the domain (e.g. around a crack tip) and use standard FE (e.g. FE_Q) elsewhere. This can be achieved by using the hp finite element framework in deal.II that allows for the use of different elements on different cells. To make the resulting space \(C^0\) continuous, it is then necessary for the DoFHandler class and DoFTools::make_hanging_node_constraints() function to be able to figure out what to do at the interface between enriched and non-enriched cells. Specifically, we want the degrees of freedom corresponding to enriched shape functions to be zero at these interfaces. These classes and functions can not to do this automatically, but the effect can be achieved by using not just a regular FE_Q on cells without enrichment, but to wrap the FE_Q into an FE_Enriched object without actually enriching it. This can be done as follows:
This constructor is equivalent to calling
and will result in the correct constraints for enriched DoFs attributed to support points on the interface between the two regions.
When using this class, please cite
The PUM was introduced in
The implementation of the class is based on FESystem which is aggregated as a private member. The simplest constructor FE_Enriched<dim> fe(FE_Q<dim>(2), FE_Q<dim>(1),function)
will internally initialize FESystem as
Note that it would not be wise to have this class derived from FESystem as the latter concatenates the given elements into different components of a vector element, whereas the current class combines the given elements into the same components. For instance, if two scalar elements are given, the resulting element will be scalar rather than have two components when doing the same with an FESystem.
The ordering of the shape function, interface_constrains
, the prolongation
(embedding) and the restriction
matrices are taken from the FESystem class.
Definition at line 216 of file fe_enriched.h.
FE_Enriched< dim, spacedim >::FE_Enriched | ( | const FiniteElement< dim, spacedim > & | fe_base, |
const FiniteElement< dim, spacedim > & | fe_enriched, | ||
const Function< spacedim > * | enrichment_function | ||
) |
Constructor which takes base FiniteElement fe_base
and the enrichment FiniteElement fe_enriched
which will be multiplied by the enrichment_function
.
In case fe_enriched
is other than FE_Nothing, the lifetime of the enrichment_function
must be at least as long as the FE_Enriched object.
Definition at line 146 of file fe_enriched.cc.
FE_Enriched< dim, spacedim >::FE_Enriched | ( | const FiniteElement< dim, spacedim > & | fe_base | ) |
Constructor which only wraps the base FE fe_base
. As for the enriched finite element space, FE_Nothing is used. Continuity constraints will be automatically generated when this non-enriched element is used in conjunction with enriched finite element within the hp::DoFHandler.
See the discussion in the class documentation on how to use this element in the context of hp finite element methods.
Definition at line 136 of file fe_enriched.cc.
FE_Enriched< dim, spacedim >::FE_Enriched | ( | const FiniteElement< dim, spacedim > * | fe_base, |
const std::vector< const FiniteElement< dim, spacedim > *> & | fe_enriched, | ||
const std::vector< std::vector< std::function< const Function< spacedim > *(const typename Triangulation< dim, spacedim >::cell_iterator &)>>> & | functions | ||
) |
Constructor which takes pointer to the base FiniteElement fe_base
and a vector of enriched FiniteElement's fe_enriched
. fe_enriched
[i] finite element will be enriched with functions in functions
[i].
This is the most general public constructor which also allows to have different enrichment functions in different disjoint parts of the domain. To that end the last argument provides an association of cell iterator to a Function. This is done to simplify the usage of this class when the number of disjoint domains with different functions is more than a few. Otherwise one would have to use different instance of this class for each disjoint enriched domain.
If you don't plan to use this feature, you can utilize C++11 lambdas to define dummy functions. Below is an example which uses two functions with the first element to be enriched and a single function with the second one.
Definition at line 165 of file fe_enriched.cc.
|
private |
The most general private constructor. The first two input parameters are consistent with those in FESystem. It is used internally only with multiplicities[0]=1
, which is a logical requirement for this finite element.
Definition at line 178 of file fe_enriched.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 284 of file fe_enriched.cc.
|
overridevirtual |
Given a set of update flags, compute which other quantities also need to be computed in order to satisfy the request by the given flags. Then return the combination of the original set of flags and those just computed.
As an example, if update_flags
contains update_gradients a finite element class will typically require the computation of the inverse of the Jacobian matrix in order to rotate the gradient of shape functions on the reference cell to the real cell. It would then return not just update_gradients, but also update_covariant_transformation, the flag that makes the mapping class produce the inverse of the Jacobian matrix.
An extensive discussion of the interaction between this function and FEValues can be found in the How Mapping, FiniteElement, and FEValues work together documentation module.
Implements FiniteElement< dim, spacedim >.
Definition at line 302 of file fe_enriched.cc.
|
overridevirtual |
Return a string that identifies a finite element.
Implements FiniteElement< dim, spacedim >.
Definition at line 486 of file fe_enriched.cc.
|
overridevirtual |
Access to a composing element. The index needs to be smaller than the number of base elements. In the context of this class, the number of base elements is always more than one: a non-enriched element plus an element to be enriched, which could be FE_Nothing.
Reimplemented from FiniteElement< dim, spacedim >.
Definition at line 507 of file fe_enriched.cc.
|
overridevirtual |
Return the value of the ith
shape function at the point p
. p
is a point on the reference element.
This function returns meaningful values only for non-enriched element as real-space enrichment requires evaluation of the function at the point in real-space.
Reimplemented from FiniteElement< dim, spacedim >.
Definition at line 271 of file fe_enriched.cc.
|
overridevirtual |
Projection from a fine grid space onto a coarse grid space.
This function only makes sense when all child elements are also enriched using the same function(s) as the parent element.
Reimplemented from FiniteElement< dim, spacedim >.
Definition at line 1025 of file fe_enriched.cc.
|
overridevirtual |
Embedding matrix between grids.
This function only makes sense when all child elements are also enriched using the same function(s) as the parent element.
Reimplemented from FiniteElement< dim, spacedim >.
Definition at line 1015 of file fe_enriched.cc.
|
overridevirtual |
Return whether this element implements hp constraints.
This function returns true
if and only if all its base elements return true
for this function.
Reimplemented from FiniteElement< dim, spacedim >.
Definition at line 873 of file fe_enriched.cc.
|
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.
Reimplemented from FiniteElement< dim, spacedim >.
Definition at line 881 of file fe_enriched.cc.
|
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.
Reimplemented from FiniteElement< dim, spacedim >.
Definition at line 903 of file fe_enriched.cc.
|
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.
Reimplemented from FiniteElement< dim, spacedim >.
Definition at line 927 of file fe_enriched.cc.
|
overridevirtual |
Same as hp_vertex_dof_indices(), except that the function treats degrees of freedom on lines.
Reimplemented from FiniteElement< dim, spacedim >.
Definition at line 945 of file fe_enriched.cc.
|
overridevirtual |
Same as hp_vertex_dof_indices(), except that the function treats degrees of freedom on quads.
Reimplemented from FiniteElement< dim, spacedim >.
Definition at line 963 of file fe_enriched.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 981 of file fe_enriched.cc.
const std::vector< std::vector< std::function< const Function< spacedim > *(const typename Triangulation< dim, spacedim >::cell_iterator &)> > > FE_Enriched< dim, spacedim >::get_enrichments | ( | ) | const |
Return enrichment functions
Definition at line 263 of file fe_enriched.cc.
const FESystem< dim, spacedim > & FE_Enriched< dim, spacedim >::get_fe_system | ( | ) | const |
Return the underlying FESystem object.
Definition at line 865 of file fe_enriched.cc.
|
protected |
Auxiliary function called from get_data, get_face_data and get_subface_data. It take internal data of FESystem object in fes_data
and the quadrature rule qudrature
.
This function essentially take the internal data from an instance of FESystem class and wraps it into our own InternalData class which additionally has objects to hold values/gradients/hessians of enrichment functions at each quadrature point depending on flags
.
Definition at line 326 of file fe_enriched.cc.
|
overrideprotectedvirtual |
Prepare internal data structures and fill in values independent of the cell. Returns a pointer to an object of which the caller of this function (FEValues) then has to assume ownership (which includes destruction when it is no more needed).
Implements FiniteElement< dim, spacedim >.
Definition at line 407 of file fe_enriched.cc.
|
overrideprotectedvirtual |
Like get_data(), but return an object that will later be used for evaluating shape function information at quadrature points on faces of cells. The object will then be used in calls to implementations of FiniteElement::fill_fe_face_values(). See the documentation of get_data() for more information.
The default implementation of this function converts the face quadrature into a cell quadrature with appropriate quadrature point locations, and with that calls the get_data() function above that has to be implemented in derived classes.
[in] | update_flags | A set of UpdateFlags values that describe what kind of information the FEValues object requests the finite element to compute. This set of flags may also include information that the finite element can not compute, e.g., flags that pertain to data produced by the mapping. An implementation of this function needs to set up all data fields in the returned object that are necessary to produce the finite- element related data specified by these flags, and may already pre- compute part of this information as discussed above. Elements may want to store these update flags (or a subset of these flags) in InternalDataBase::update_each so they know at the time when FiniteElement::fill_fe_face_values() is called what they are supposed to compute |
[in] | mapping | A reference to the mapping used for computing values and derivatives of shape functions. |
[in] | quadrature | A reference to the object that describes where the shape functions should be evaluated. |
[out] | output_data | A reference to the object that FEValues will use in conjunction with the object returned here and where an implementation of FiniteElement::fill_fe_face_values() will place the requested information. This allows the current function to already pre-compute pieces of information that can be computed on the reference cell, as discussed above. FEValues guarantees that this output object and the object returned by the current function will always be used together. |
Reimplemented from FiniteElement< dim, spacedim >.
Definition at line 368 of file fe_enriched.cc.
|
overrideprotectedvirtual |
Like get_data(), but return an object that will later be used for evaluating shape function information at quadrature points on children of faces of cells. The object will then be used in calls to implementations of FiniteElement::fill_fe_subface_values(). See the documentation of get_data() for more information.
The default implementation of this function converts the face quadrature into a cell quadrature with appropriate quadrature point locations, and with that calls the get_data() function above that has to be implemented in derived classes.
[in] | update_flags | A set of UpdateFlags values that describe what kind of information the FEValues object requests the finite element to compute. This set of flags may also include information that the finite element can not compute, e.g., flags that pertain to data produced by the mapping. An implementation of this function needs to set up all data fields in the returned object that are necessary to produce the finite- element related data specified by these flags, and may already pre- compute part of this information as discussed above. Elements may want to store these update flags (or a subset of these flags) in InternalDataBase::update_each so they know at the time when FiniteElement::fill_fe_subface_values() is called what they are supposed to compute |
[in] | mapping | A reference to the mapping used for computing values and derivatives of shape functions. |
[in] | quadrature | A reference to the object that describes where the shape functions should be evaluated. |
[out] | output_data | A reference to the object that FEValues will use in conjunction with the object returned here and where an implementation of FiniteElement::fill_fe_subface_values() will place the requested information. This allows the current function to already pre-compute pieces of information that can be computed on the reference cell, as discussed above. FEValues guarantees that this output object and the object returned by the current function will always be used together. |
Reimplemented from FiniteElement< dim, spacedim >.
Definition at line 387 of file fe_enriched.cc.
|
private |
This function sets up the index table for the system as well as restriction
and prolongation
matrices.
Definition at line 425 of file fe_enriched.cc.
|
private |
After calling fill_fe_(face/subface_)values this function implements the chain rule to multiply stored shape values/gradient/hessians by those of enrichment function evaluated at quadrature points.
Definition at line 623 of file fe_enriched.cc.
|
protected |
For each finite element i
used in enrichment and each enrichment function j
associated with it (essentially its multiplicity), base_no_mult_local_enriched_dofs
[i][j] contains the associated local DoFs on the FE_Enriched finite element.
Definition at line 559 of file fe_enriched.h.
|
protected |
Enrichment functions. The size of the first vector is the same as the number of FiniteElement spaces used with enrichment. Whereas the size of the inner vector corresponds to the number of enrichment functions associated with a single FiniteElement.
Definition at line 570 of file fe_enriched.h.
|
protected |
Auxiliary variable used to distinguish between the case when we do enrichment and when the class simply wraps another FiniteElement.
This variable is initialized in the constructor by looping over a vector of enrichment elements and checking if all of them are FE_Nothing. If this is the case, then the value is set to false
, otherwise it is true
.
Definition at line 581 of file fe_enriched.h.
|
private |
The underlying FESystem object.
Definition at line 694 of file fe_enriched.h.