16 #ifndef dealii_fe_nedelec_sz_h 17 #define dealii_fe_nedelec_sz_h 19 #include <deal.II/base/derivative_form.h> 20 #include <deal.II/base/polynomials_integrated_legendre_sz.h> 21 #include <deal.II/base/qprojector.h> 23 #include <deal.II/fe/fe.h> 24 #include <deal.II/fe/fe_values.h> 25 #include <deal.II/fe/mapping.h> 27 DEAL_II_NAMESPACE_OPEN
72 template <
int dim,
int spacedim = dim>
76 static_assert(dim == spacedim,
77 "FE_NedelecSZ is only implemented for dim==spacedim!");
90 virtual std::unique_ptr<FiniteElement<dim, dim>>
91 clone()
const override;
106 const unsigned int component)
const override;
121 const unsigned int component)
const override;
136 const unsigned int component)
const override;
161 virtual std::unique_ptr<
162 typename ::FiniteElement<dim, spacedim>::InternalDataBase>
169 &output_data)
const override;
183 const ::internal::FEValuesImplementation::MappingRelatedData<dim, dim>
187 &data)
const override;
197 const unsigned int face_no,
201 const ::internal::FEValuesImplementation::MappingRelatedData<dim, dim>
205 &data)
const override;
213 const unsigned int face_no,
214 const unsigned int sub_no,
218 const ::internal::FEValuesImplementation::MappingRelatedData<dim, dim>
222 &data)
const override;
261 mutable std::vector<std::vector<DerivativeForm<1, dim, dim>>>
shape_grads;
413 static std::vector<unsigned int>
457 DEAL_II_NAMESPACE_CLOSE
virtual void fill_fe_values(const typename Triangulation< dim, dim >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const Quadrature< dim > &quadrature, const Mapping< dim, dim > &mapping, const typename Mapping< dim, dim >::InternalDataBase &mapping_internal, const ::internal::FEValuesImplementation::MappingRelatedData< dim, dim > &mapping_data, const typename FiniteElement< dim, dim >::InternalDataBase &fedata, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, dim > &data) const override
std::vector< std::vector< double > > edge_lambda_grads_2d
virtual double shape_value_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
UpdateFlags update_once(const UpdateFlags flags) const
std::vector< std::vector< double > > face_lambda_values
void fill_edge_values(const typename Triangulation< dim, dim >::cell_iterator &cell, const Quadrature< dim > &quadrature, const InternalData &fedata) const
std::vector< std::vector< double > > edge_sigma_values
const unsigned int degree
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
void fill_face_values(const typename Triangulation< dim, dim >::cell_iterator &cell, const Quadrature< dim > &quadrature, const InternalData &fedata) const
virtual double shape_value(const unsigned int i, const Point< dim > &p) const override
virtual Tensor< 1, dim > shape_grad(const unsigned int i, const Point< dim > &p) const override
virtual std::string get_name() const override
virtual void fill_fe_face_values(const typename Triangulation< dim, dim >::cell_iterator &cell, const unsigned int face_no, const Quadrature< dim - 1 > &quadrature, const Mapping< dim, dim > &mapping, const typename Mapping< dim, dim >::InternalDataBase &mapping_internal, const ::internal::FEValuesImplementation::MappingRelatedData< dim, dim > &mapping_data, const typename FiniteElement< dim, dim >::InternalDataBase &fedata, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, dim > &data) const override
std::vector< std::vector< DerivativeForm< 1, dim, dim > > > shape_grads
std::vector< std::vector< std::vector< double > > > edge_lambda_gradgrads_3d
virtual void fill_fe_subface_values(const typename Triangulation< dim, dim >::cell_iterator &cell, const unsigned int face_no, const unsigned int sub_no, const Quadrature< dim - 1 > &quadrature, const Mapping< dim, dim > &mapping, const typename Mapping< dim, dim >::InternalDataBase &mapping_internal, const ::internal::FEValuesImplementation::MappingRelatedData< dim, dim > &mapping_data, const typename FiniteElement< dim, dim >::InternalDataBase &fedata, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, dim > &data) const override
Abstract base class for mapping classes.
std::vector< Polynomials::Polynomial< double > > IntegratedLegendrePolynomials
std::vector< std::vector< Tensor< 1, dim > > > shape_values
virtual Tensor< 2, dim > shape_grad_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 std::unique_ptr< FiniteElement< dim, dim > > clone() const override
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
void create_polynomials(const unsigned int degree)
std::vector< std::vector< double > > edge_lambda_values
UpdateFlags update_each(const UpdateFlags flags) const
virtual Tensor< 1, dim > shape_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
std::vector< std::vector< std::vector< double > > > edge_lambda_grads_3d
FE_NedelecSZ(const unsigned int degree)
std::vector< std::vector< double > > edge_sigma_grads
std::vector< std::vector< std::vector< double > > > sigma_imj_values
std::vector< std::vector< double > > face_lambda_grads
virtual std::unique_ptr< typename ::FiniteElement< dim, spacedim >::InternalDataBase > get_data(const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
std::vector< std::vector< std::vector< double > > > sigma_imj_grads
unsigned int compute_num_dofs(const unsigned int degree) const