 |
Reference documentation for deal.II version 9.1.1
|
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
16 #ifndef dealii_fe_p1nc_h
17 #define dealii_fe_p1nc_h
19 #include <deal.II/base/config.h>
21 #include <deal.II/base/qprojector.h>
22 #include <deal.II/base/quadrature.h>
24 #include <deal.II/fe/fe.h>
27 DEAL_II_NAMESPACE_OPEN
287 virtual std::unique_ptr<FiniteElement<2, 2>>
288 clone()
const override;
293 virtual ~FE_P1NC()
override =
default;
302 static std::vector<unsigned int>
310 static std::array<std::array<double, 3>, 4>
321 virtual std::unique_ptr<FiniteElement<2, 2>::InternalDataBase>
323 const UpdateFlags update_flags,
327 &output_data)
const override;
329 virtual std::unique_ptr<FiniteElement<2, 2>::InternalDataBase>
331 const UpdateFlags update_flags,
335 &output_data)
const override;
337 virtual std::unique_ptr<FiniteElement<2, 2>::InternalDataBase>
339 const UpdateFlags update_flags,
343 &output_data)
const override;
359 &output_data)
const override;
367 const unsigned int face_no,
371 const ::internal::FEValuesImplementation::MappingRelatedData<2, 2>
373 const InternalDataBase &fe_internal,
375 &output_data)
const override;
383 const unsigned int face_no,
384 const unsigned int sub_no,
388 const ::internal::FEValuesImplementation::MappingRelatedData<2, 2>
390 const InternalDataBase &fe_internal,
392 &output_data)
const override;
406 DEAL_II_NAMESPACE_CLOSE
virtual void fill_fe_values(const Triangulation< 2, 2 >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const Quadrature< 2 > &quadrature, const Mapping< 2, 2 > &mapping, const Mapping< 2, 2 >::InternalDataBase &mapping_internal, const internal::FEValuesImplementation::MappingRelatedData< 2, 2 > &mapping_data, const FiniteElement< 2, 2 >::InternalDataBase &fe_internal, internal::FEValuesImplementation::FiniteElementRelatedData< 2, 2 > &output_data) const override
virtual std::unique_ptr< FiniteElement< 2, 2 > > clone() const override
virtual UpdateFlags requires_update_flags(const UpdateFlags flags) const override
virtual ~FE_P1NC() override=default
Abstract base class for mapping classes.
static std::vector< unsigned int > get_dpo_vector()
virtual std::string get_name() const override
void initialize_constraints()
static std::array< std::array< double, 3 >, 4 > get_linear_shape_coefficients(const Triangulation< 2, 2 >::cell_iterator &cell)
virtual void fill_fe_face_values(const Triangulation< 2, 2 >::cell_iterator &cell, const unsigned int face_no, const Quadrature< 1 > &quadrature, const Mapping< 2, 2 > &mapping, const Mapping< 2, 2 >::InternalDataBase &mapping_internal, const ::internal::FEValuesImplementation::MappingRelatedData< 2, 2 > &mapping_data, const InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< 2, 2 > &output_data) const override
virtual void fill_fe_subface_values(const Triangulation< 2, 2 >::cell_iterator &cell, const unsigned int face_no, const unsigned int sub_no, const Quadrature< 1 > &quadrature, const Mapping< 2, 2 > &mapping, const Mapping< 2, 2 >::InternalDataBase &mapping_internal, const ::internal::FEValuesImplementation::MappingRelatedData< 2, 2 > &mapping_data, const InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< 2, 2 > &output_data) const override
virtual std::unique_ptr< FiniteElement< 2, 2 >::InternalDataBase > get_data(const UpdateFlags update_flags, const Mapping< 2, 2 > &, const Quadrature< 2 > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< 2, 2 > &output_data) const override