16#ifndef dealii_fe_face_h
17#define dealii_fe_face_h
54template <
int dim,
int spacedim = dim>
56 :
public FE_PolyFace<TensorProductPolynomials<dim - 1>, dim, spacedim>
66 virtual std::unique_ptr<FiniteElement<dim, spacedim>>
67 clone()
const override;
87 std::vector<double> & nodal_values)
const override;
100 const unsigned int face_no = 0)
const override;
113 const unsigned int subface,
115 const unsigned int face_no = 0)
const override;
123 const unsigned int face_index)
const override;
147 virtual std::vector<std::pair<unsigned int, unsigned int>>
157 virtual std::vector<std::pair<unsigned int, unsigned int>>
167 virtual std::vector<std::pair<unsigned int, unsigned int>>
169 const unsigned int face_no = 0)
const override;
183 const unsigned int codim = 0)
const override final;
193 virtual std::pair<Table<2, bool>, std::vector<unsigned int>>
200 static std::vector<unsigned int>
216template <
int spacedim>
225 virtual std::unique_ptr<FiniteElement<1, spacedim>>
226 clone()
const override;
251 const unsigned int face_no = 0)
const override;
264 const unsigned int subface,
266 const unsigned int face_no = 0)
const override;
274 const unsigned int face_index)
const override;
300 virtual std::vector<std::pair<unsigned int, unsigned int>>
310 virtual std::vector<std::pair<unsigned int, unsigned int>>
320 virtual std::vector<std::pair<unsigned int, unsigned int>>
322 const unsigned int face_no = 0)
const override;
328 virtual std::pair<Table<2, bool>, std::vector<unsigned int>>
339 virtual std::unique_ptr<typename FiniteElement<1, spacedim>::InternalDataBase>
348 return std::make_unique<
354 std::unique_ptr<typename FiniteElement<1, spacedim>::InternalDataBase>
367 std::make_unique<typename FiniteElement<1, spacedim>::InternalDataBase>();
370 const unsigned int n_q_points = quadrature[0].
size();
384 std::unique_ptr<typename FiniteElement<1, spacedim>::InternalDataBase>
391 &output_data)
const override
411 &output_data)
const override;
418 const unsigned int face_no,
427 &output_data)
const override;
432 const unsigned int face_no,
433 const unsigned int sub_no,
442 &output_data)
const override;
448 static std::vector<unsigned int>
476template <
int dim,
int spacedim = dim>
487 virtual std::unique_ptr<FiniteElement<dim, spacedim>>
488 clone()
const override;
509 const unsigned int face_no = 0)
const override;
522 const unsigned int subface,
524 const unsigned int face_no = 0)
const override;
532 const unsigned int face_index)
const override;
546 const unsigned int codim = 0)
const override final;
554 virtual std::pair<Table<2, bool>, std::vector<unsigned int>>
561 static std::vector<unsigned int>
571template <
int spacedim>
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const override
virtual bool hp_constraints_are_implemented() const override
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const override
virtual FiniteElementDomination::Domination compare_for_domination(const FiniteElement< dim, spacedim > &fe_other, const unsigned int codim=0) const override final
virtual std::string get_name() const override
static std::vector< unsigned int > get_dpo_vector(const unsigned int deg)
virtual void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix, const unsigned int face_no=0) const override
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const override
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix, const unsigned int face_no=0) const override
std::unique_ptr< typename FiniteElement< 1, spacedim >::InternalDataBase > get_face_data(const UpdateFlags update_flags, const Mapping< 1, spacedim > &, const hp::QCollection< 0 > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< 1, spacedim > &) const override
virtual std::unique_ptr< typename FiniteElement< 1, spacedim >::InternalDataBase > get_data(const UpdateFlags, const Mapping< 1, spacedim > &, const Quadrature< 1 > &, ::internal::FEValuesImplementation::FiniteElementRelatedData< 1, spacedim > &) const override
std::unique_ptr< typename FiniteElement< 1, spacedim >::InternalDataBase > get_subface_data(const UpdateFlags update_flags, const Mapping< 1, spacedim > &mapping, const Quadrature< 0 > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< 1, spacedim > &output_data) const override
virtual FiniteElementDomination::Domination compare_for_domination(const FiniteElement< dim, spacedim > &fe_other, const unsigned int codim=0) const override final
virtual std::string get_name() const override
static std::vector< unsigned int > get_dpo_vector(const unsigned int deg)
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::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const override
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix, const unsigned int face_no=0) 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::unique_ptr< FiniteElement< dim, spacedim > > clone() const override
virtual void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix, const unsigned int face_no=0) const override
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) 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_quad_dof_identities(const FiniteElement< dim, spacedim > &fe_other, const unsigned int face_no=0) const override
virtual bool hp_constraints_are_implemented() const override
virtual void fill_fe_face_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const hp::QCollection< dim - 1 > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const typename FiniteElement< dim, spacedim >::InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
std::unique_ptr< typename FiniteElement< dim, spacedim >::InternalDataBase > get_face_data(const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const hp::QCollection< dim - 1 > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
virtual void fill_fe_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const Quadrature< dim > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const typename FiniteElement< dim, spacedim >::InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
virtual void fill_fe_subface_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int sub_no, const Quadrature< dim - 1 > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const typename FiniteElement< dim, spacedim >::InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
Abstract base class for mapping classes.
unsigned int size() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
@ update_hessians
Second derivatives of shape functions.
@ update_gradients
Shape function gradients.