20#include <deal.II/fe/fe_poly_face.templates.h>
33 namespace FE_FaceQImplementation
38 get_QGaussLobatto_points(
const unsigned int degree)
43 return std::vector<Point<1>>(1,
Point<1>(0.5));
49template <
int dim,
int spacedim>
54 internal::FE_FaceQImplementation::get_QGaussLobatto_points(degree))),
62 const unsigned int codim = dim - 1;
64 Utilities::fixed_power<codim>(this->degree + 1));
66 if (this->degree == 0)
67 for (
unsigned int d = 0; d < codim; ++d)
71 std::vector<Point<1>> points =
72 internal::FE_FaceQImplementation::get_QGaussLobatto_points(
degree);
75 for (
unsigned int iz = 0; iz <= ((codim > 2) ? this->degree : 0); ++iz)
76 for (
unsigned int iy = 0; iy <= ((codim > 1) ? this->degree : 0); ++iy)
77 for (
unsigned int ix = 0; ix <= this->
degree; ++ix)
97 for (
unsigned int i = 0; i < n_face_dofs; ++i)
98 for (
unsigned int d = 0; d < dim; ++d)
100 for (
unsigned int e = 0, c = 0; e < dim; ++e)
104 unsigned int renumber = i;
105 if (dim == 3 && d == 1)
119template <
int dim,
int spacedim>
120std::unique_ptr<FiniteElement<dim, spacedim>>
123 return std::make_unique<FE_FaceQ<dim, spacedim>>(this->degree);
128template <
int dim,
int spacedim>
135 std::ostringstream namebuf;
137 << this->degree <<
")";
139 return namebuf.str();
144template <
int dim,
int spacedim>
149 const unsigned int face_no)
const
151 get_subface_interpolation_matrix(source_fe,
153 interpolation_matrix,
159template <
int dim,
int spacedim>
163 const unsigned int subface,
165 const unsigned int face_no)
const
169 Assert(interpolation_matrix.
n() == this->n_dofs_per_face(face_no),
171 this->n_dofs_per_face(face_no)));
186 this->n_dofs_per_face(face_no) <= source_fe->n_dofs_per_face(face_no),
188 spacedim>::ExcInterpolationNotImplemented()));
192 source_fe->get_unit_face_support_points(face_no));
197 const double eps = 2e-13 * (this->degree + 1) * (dim - 1);
201 for (
unsigned int i = 0; i < source_fe->n_dofs_per_face(face_no); ++i)
203 const Point<dim - 1> p =
205 face_quadrature.point(i) :
207 face_quadrature.point(i), subface);
209 for (
unsigned int j = 0; j < this->n_dofs_per_face(face_no); ++j)
211 double matrix_entry = this->poly_space.compute_value(j, p);
216 if (std::fabs(matrix_entry - 1.0) < eps)
218 if (std::fabs(matrix_entry) < eps)
221 interpolation_matrix(i, j) = matrix_entry;
228 for (
unsigned int j = 0; j < source_fe->n_dofs_per_face(face_no); ++j)
232 for (
unsigned int i = 0; i < this->n_dofs_per_face(face_no); ++i)
233 sum += interpolation_matrix(j, i);
239 else if (
dynamic_cast<const FE_Nothing<dim> *
>(&x_source_fe) !=
nullptr)
247 spacedim>::ExcInterpolationNotImplemented()));
252template <
int dim,
int spacedim>
255 const unsigned int shape_index,
256 const unsigned int face_index)
const
258 return (face_index == (shape_index / this->n_dofs_per_face(face_index)));
263template <
int dim,
int spacedim>
264std::vector<unsigned int>
267 std::vector<unsigned int> dpo(dim + 1, 0U);
268 dpo[dim - 1] = deg + 1;
269 for (
unsigned int i = 1; i < dim - 1; ++i)
270 dpo[dim - 1] *= deg + 1;
276template <
int dim,
int spacedim>
285template <
int dim,
int spacedim>
286std::vector<std::pair<unsigned int, unsigned int>>
291 return std::vector<std::pair<unsigned int, unsigned int>>();
296template <
int dim,
int spacedim>
297std::vector<std::pair<unsigned int, unsigned int>>
305 return std::vector<std::pair<unsigned int, unsigned int>>();
318 const unsigned int p = this->degree;
319 const unsigned int q = fe_q_other->degree;
321 std::vector<std::pair<unsigned int, unsigned int>> identities;
323 const std::vector<unsigned int> &index_map_inverse =
324 this->poly_space.get_numbering_inverse();
325 const std::vector<unsigned int> &index_map_inverse_other =
326 fe_q_other->poly_space.get_numbering_inverse();
328 for (
unsigned int i = 0; i < p + 1; ++i)
329 for (
unsigned int j = 0; j < q + 1; ++j)
331 this->unit_support_points[index_map_inverse[i]][dim - 1] -
332 fe_q_other->unit_support_points[index_map_inverse_other[j]]
334 identities.emplace_back(i, j);
342 return std::vector<std::pair<unsigned int, unsigned int>>();
354 return std::vector<std::pair<unsigned int, unsigned int>>();
359 return std::vector<std::pair<unsigned int, unsigned int>>();
366template <
int dim,
int spacedim>
367std::vector<std::pair<unsigned int, unsigned int>>
370 const unsigned int)
const
376 return std::vector<std::pair<unsigned int, unsigned int>>();
388 const unsigned int p = this->degree;
389 const unsigned int q = fe_q_other->degree;
391 std::vector<std::pair<unsigned int, unsigned int>> identities;
393 const std::vector<unsigned int> &index_map_inverse =
394 this->poly_space.get_numbering_inverse();
395 const std::vector<unsigned int> &index_map_inverse_other =
396 fe_q_other->poly_space.get_numbering_inverse();
398 std::vector<std::pair<unsigned int, unsigned int>> identities_1d;
400 for (
unsigned int i = 0; i < p + 1; ++i)
401 for (
unsigned int j = 0; j < q + 1; ++j)
403 this->unit_support_points[index_map_inverse[i]][dim - 2] -
404 fe_q_other->unit_support_points[index_map_inverse_other[j]]
406 identities_1d.emplace_back(i, j);
408 for (
unsigned int n1 = 0; n1 < identities_1d.size(); ++n1)
409 for (
unsigned int n2 = 0; n2 < identities_1d.size(); ++n2)
410 identities.emplace_back(identities_1d[n1].first * (p + 1) +
411 identities_1d[n2].first,
412 identities_1d[n1].second * (q + 1) +
413 identities_1d[n2].
second);
421 return std::vector<std::pair<unsigned int, unsigned int>>();
433 return std::vector<std::pair<unsigned int, unsigned int>>();
438 return std::vector<std::pair<unsigned int, unsigned int>>();
445template <
int dim,
int spacedim>
449 const unsigned int codim)
const
458 if (this->degree < fe_faceq_other->degree)
460 else if (this->degree == fe_faceq_other->degree)
468 if (fe_nothing->is_dominating())
481template <
int dim,
int spacedim>
482std::pair<Table<2, bool>, std::vector<unsigned int>>
486 for (
unsigned int i = 0; i < this->n_dofs_per_cell(); ++i)
487 constant_modes(0, i) =
true;
488 return std::pair<Table<2, bool>, std::vector<unsigned int>>(
489 constant_modes, std::vector<unsigned int>(1, 0));
492template <
int dim,
int spacedim>
496 std::vector<double> &nodal_values)
const
499 this->get_unit_support_points().size());
503 for (
unsigned int i = 0; i < this->n_dofs_per_cell(); ++i)
507 nodal_values[i] = support_point_values[i](0);
513template <
int spacedim>
533template <
int spacedim>
534std::unique_ptr<FiniteElement<1, spacedim>>
537 return std::make_unique<FE_FaceQ<1, spacedim>>(this->degree);
542template <
int spacedim>
549 std::ostringstream namebuf;
551 << this->degree <<
")";
553 return namebuf.str();
558template <
int spacedim>
563 const unsigned int face_no)
const
565 get_subface_interpolation_matrix(source_fe,
567 interpolation_matrix,
573template <
int spacedim>
579 const unsigned int face_no)
const
581 Assert(interpolation_matrix.
n() == this->n_dofs_per_face(face_no),
583 this->n_dofs_per_face(face_no)));
587 interpolation_matrix(0, 0) = 1.;
592template <
int spacedim>
595 const unsigned int face_index)
const
598 return (face_index == shape_index);
603template <
int spacedim>
604std::vector<unsigned int>
607 std::vector<unsigned int> dpo(2, 0U);
614template <
int spacedim>
621template <
int spacedim>
622std::vector<std::pair<unsigned int, unsigned int>>
627 return std::vector<std::pair<unsigned int, unsigned int>>(1,
634template <
int spacedim>
635std::vector<std::pair<unsigned int, unsigned int>>
640 return std::vector<std::pair<unsigned int, unsigned int>>();
645template <
int spacedim>
646std::vector<std::pair<unsigned int, unsigned int>>
649 const unsigned int)
const
652 return std::vector<std::pair<unsigned int, unsigned int>>();
657template <
int spacedim>
658std::pair<Table<2, bool>, std::vector<unsigned int>>
662 for (
unsigned int i = 0; i < this->n_dofs_per_cell(); ++i)
663 constant_modes(0, i) =
true;
664 return std::pair<Table<2, bool>, std::vector<unsigned int>>(
665 constant_modes, std::vector<unsigned int>(1, 0));
670template <
int spacedim>
686template <
int spacedim>
705template <
int spacedim>
709 const unsigned int face,
719 const unsigned int foffset = face;
722 for (
unsigned int k = 0; k < this->n_dofs_per_cell(); ++k)
723 output_data.shape_values(k, 0) = 0.;
724 output_data.shape_values(foffset, 0) = 1;
729template <
int spacedim>
751template <
int dim,
int spacedim>
755 Polynomials::Legendre::generate_complete_basis(degree)),
765template <
int dim,
int spacedim>
766std::unique_ptr<FiniteElement<dim, spacedim>>
769 return std::make_unique<FE_FaceP<dim, spacedim>>(this->degree);
774template <
int dim,
int spacedim>
781 std::ostringstream namebuf;
783 << this->degree <<
")";
785 return namebuf.str();
790template <
int dim,
int spacedim>
793 const unsigned int shape_index,
794 const unsigned int face_index)
const
796 return (face_index == (shape_index / this->n_dofs_per_face(face_index)));
801template <
int dim,
int spacedim>
802std::vector<unsigned int>
805 std::vector<unsigned int> dpo(dim + 1, 0U);
806 dpo[dim - 1] = deg + 1;
807 for (
unsigned int i = 1; i < dim - 1; ++i)
809 dpo[dim - 1] *= deg + 1 + i;
810 dpo[dim - 1] /= i + 1;
817template <
int dim,
int spacedim>
826template <
int dim,
int spacedim>
830 const unsigned int codim)
const
839 if (this->degree < fe_facep_other->degree)
841 else if (this->degree == fe_facep_other->degree)
849 if (fe_nothing->is_dominating())
864template <
int dim,
int spacedim>
869 const unsigned int face_no)
const
871 get_subface_interpolation_matrix(source_fe,
873 interpolation_matrix,
879template <
int dim,
int spacedim>
883 const unsigned int subface,
885 const unsigned int face_no)
const
889 Assert(interpolation_matrix.
n() == this->n_dofs_per_face(face_no),
891 this->n_dofs_per_face(face_no)));
906 this->n_dofs_per_face(face_no) <= source_fe->n_dofs_per_face(face_no),
908 spacedim>::ExcInterpolationNotImplemented()));
913 const QGauss<dim - 1> face_quadrature(source_fe->degree + 1);
918 const double eps = 2e-13 * (this->degree + 1) * (dim - 1);
921 source_fe->n_dofs_per_face(face_no));
923 for (
unsigned int k = 0; k < face_quadrature.size(); ++k)
925 const Point<dim - 1> p =
927 face_quadrature.point(k) :
929 face_quadrature.point(k), subface);
931 for (
unsigned int j = 0; j < source_fe->n_dofs_per_face(face_no); ++j)
932 mass(k, j) = source_fe->poly_space.compute_value(j, p);
942 for (
unsigned int i = 0; i < this->n_dofs_per_face(face_no); ++i)
944 for (
unsigned int k = 0; k < face_quadrature.size(); ++k)
946 const Point<dim - 1> p =
948 face_quadrature.point(k) :
950 face_quadrature.point(k), subface);
951 v_in(k) = this->poly_space.compute_value(i, p);
953 [[maybe_unused]]
const double result = H.
least_squares(v_out, v_in);
956 for (
unsigned int j = 0; j < source_fe->n_dofs_per_face(face_no); ++j)
958 double matrix_entry = v_out(j);
963 if (std::fabs(matrix_entry - 1.0) < eps)
965 if (std::fabs(matrix_entry) < eps)
968 interpolation_matrix(j, i) = matrix_entry;
972 else if (
dynamic_cast<const FE_Nothing<dim> *
>(&x_source_fe) !=
nullptr)
980 spacedim>::ExcInterpolationNotImplemented()));
985template <
int dim,
int spacedim>
986std::pair<Table<2, bool>, std::vector<unsigned int>>
991 constant_modes(0, face * this->n_dofs_per_face(face)) =
true;
992 return std::pair<Table<2, bool>, std::vector<unsigned int>>(
993 constant_modes, std::vector<unsigned int>(1, 0));
998template <
int spacedim>
1005template <
int spacedim>
1012 std::ostringstream namebuf;
1014 << this->degree <<
")";
1016 return namebuf.str();
1022#include "fe_face.inst"
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
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
FE_FaceQ(const unsigned int p)
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
const unsigned int degree
unsigned int n_dofs_per_face(unsigned int face_no=0, unsigned int child=0) const
unsigned int n_unique_faces() 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
std::vector< std::vector< Point< dim - 1 > > > unit_face_support_points
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 InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const
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
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const =0
std::vector< Point< dim > > unit_support_points
double least_squares(Vector< number2 > &dst, const Vector< number2 > &src) const
Abstract base class for mapping classes.
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
static ::ExceptionBase & ExcLeastSquaresError(double arg1)
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
@ update_hessians
Second derivatives of shape functions.
@ update_values
Shape function values.
@ update_normal_vectors
Normal vectors.
@ update_JxW_values
Transformed quadrature weights.
@ update_covariant_transformation
Covariant transformation.
@ update_gradients
Shape function gradients.
#define DEAL_II_NOT_IMPLEMENTED()
@ either_element_can_dominate
@ other_element_dominates
@ neither_element_dominates
std::string dim_string(const int dim, const int spacedim)
static const unsigned int invalid_unsigned_int
static std_cxx20::ranges::iota_view< unsigned int, unsigned int > face_indices()
static Point< dim > child_to_cell_coordinates(const Point< dim > &p, const unsigned int child_index, const RefinementCase< dim > refine_case=RefinementCase< dim >::isotropic_refinement)