21#include <deal.II/fe/fe_poly_face.templates.h>
34 namespace FE_FaceQImplementation
39 get_QGaussLobatto_points(
const unsigned int degree)
44 return std::vector<Point<1>>(1,
Point<1>(0.5));
50template <
int dim,
int spacedim>
55 internal::FE_FaceQImplementation::get_QGaussLobatto_points(degree))),
63 const unsigned int codim = dim - 1;
65 Utilities::fixed_power<codim>(this->degree + 1));
67 if (this->degree == 0)
68 for (
unsigned int d = 0; d < codim; ++d)
72 std::vector<Point<1>> points =
73 internal::FE_FaceQImplementation::get_QGaussLobatto_points(
degree);
76 for (
unsigned int iz = 0; iz <= ((codim > 2) ? this->degree : 0); ++iz)
77 for (
unsigned int iy = 0; iy <= ((codim > 1) ? this->degree : 0); ++iy)
78 for (
unsigned int ix = 0; ix <= this->
degree; ++ix)
98 for (
unsigned int i = 0; i < n_face_dofs; ++i)
99 for (
unsigned int d = 0; d < dim; ++d)
101 for (
unsigned int e = 0, c = 0; e < dim; ++e)
105 unsigned int renumber = i;
106 if (dim == 3 && d == 1)
120template <
int dim,
int spacedim>
121std::unique_ptr<FiniteElement<dim, spacedim>>
124 return std::make_unique<FE_FaceQ<dim, spacedim>>(this->degree);
129template <
int dim,
int spacedim>
136 std::ostringstream namebuf;
138 << this->degree <<
")";
140 return namebuf.str();
145template <
int dim,
int spacedim>
150 const unsigned int face_no)
const
152 get_subface_interpolation_matrix(source_fe,
154 interpolation_matrix,
160template <
int dim,
int spacedim>
164 const unsigned int subface,
166 const unsigned int face_no)
const
170 Assert(interpolation_matrix.
n() == this->n_dofs_per_face(face_no),
172 this->n_dofs_per_face(face_no)));
187 this->n_dofs_per_face(face_no) <= source_fe->n_dofs_per_face(face_no),
189 spacedim>::ExcInterpolationNotImplemented()));
193 source_fe->get_unit_face_support_points(face_no));
198 const double eps = 2e-13 * (this->degree + 1) * (dim - 1);
202 for (
unsigned int i = 0; i < source_fe->n_dofs_per_face(face_no); ++i)
204 const Point<dim - 1> p =
206 face_quadrature.point(i) :
208 face_quadrature.point(i), subface);
210 for (
unsigned int j = 0; j < this->n_dofs_per_face(face_no); ++j)
212 double matrix_entry = this->poly_space.compute_value(j, p);
217 if (std::fabs(matrix_entry - 1.0) < eps)
219 if (std::fabs(matrix_entry) < eps)
222 interpolation_matrix(i, j) = matrix_entry;
229 for (
unsigned int j = 0; j < source_fe->n_dofs_per_face(face_no); ++j)
233 for (
unsigned int i = 0; i < this->n_dofs_per_face(face_no); ++i)
234 sum += interpolation_matrix(j, i);
240 else if (
dynamic_cast<const FE_Nothing<dim> *
>(&x_source_fe) !=
nullptr)
248 spacedim>::ExcInterpolationNotImplemented()));
253template <
int dim,
int spacedim>
256 const unsigned int shape_index,
257 const unsigned int face_index)
const
259 return (face_index == (shape_index / this->n_dofs_per_face(face_index)));
264template <
int dim,
int spacedim>
265std::vector<unsigned int>
268 std::vector<unsigned int> dpo(dim + 1, 0U);
269 dpo[dim - 1] = deg + 1;
270 for (
unsigned int i = 1; i < dim - 1; ++i)
271 dpo[dim - 1] *= deg + 1;
277template <
int dim,
int spacedim>
286template <
int dim,
int spacedim>
287std::vector<std::pair<unsigned int, unsigned int>>
292 return std::vector<std::pair<unsigned int, unsigned int>>();
297template <
int dim,
int spacedim>
298std::vector<std::pair<unsigned int, unsigned int>>
306 return std::vector<std::pair<unsigned int, unsigned int>>();
319 const unsigned int p = this->degree;
320 const unsigned int q = fe_q_other->degree;
322 std::vector<std::pair<unsigned int, unsigned int>> identities;
324 const std::vector<unsigned int> &index_map_inverse =
325 this->poly_space.get_numbering_inverse();
326 const std::vector<unsigned int> &index_map_inverse_other =
327 fe_q_other->poly_space.get_numbering_inverse();
329 for (
unsigned int i = 0; i < p + 1; ++i)
330 for (
unsigned int j = 0; j < q + 1; ++j)
332 this->unit_support_points[index_map_inverse[i]][dim - 1] -
333 fe_q_other->unit_support_points[index_map_inverse_other[j]]
335 identities.emplace_back(i, j);
343 return std::vector<std::pair<unsigned int, unsigned int>>();
355 return std::vector<std::pair<unsigned int, unsigned int>>();
360 return std::vector<std::pair<unsigned int, unsigned int>>();
367template <
int dim,
int spacedim>
368std::vector<std::pair<unsigned int, unsigned int>>
371 const unsigned int)
const
377 return std::vector<std::pair<unsigned int, unsigned int>>();
389 const unsigned int p = this->degree;
390 const unsigned int q = fe_q_other->degree;
392 std::vector<std::pair<unsigned int, unsigned int>> identities;
394 const std::vector<unsigned int> &index_map_inverse =
395 this->poly_space.get_numbering_inverse();
396 const std::vector<unsigned int> &index_map_inverse_other =
397 fe_q_other->poly_space.get_numbering_inverse();
399 std::vector<std::pair<unsigned int, unsigned int>> identities_1d;
401 for (
unsigned int i = 0; i < p + 1; ++i)
402 for (
unsigned int j = 0; j < q + 1; ++j)
404 this->unit_support_points[index_map_inverse[i]][dim - 2] -
405 fe_q_other->unit_support_points[index_map_inverse_other[j]]
407 identities_1d.emplace_back(i, j);
409 for (
unsigned int n1 = 0; n1 < identities_1d.size(); ++n1)
410 for (
unsigned int n2 = 0; n2 < identities_1d.size(); ++n2)
411 identities.emplace_back(identities_1d[n1].first * (p + 1) +
412 identities_1d[n2].first,
413 identities_1d[n1].second * (q + 1) +
414 identities_1d[n2].
second);
422 return std::vector<std::pair<unsigned int, unsigned int>>();
434 return std::vector<std::pair<unsigned int, unsigned int>>();
439 return std::vector<std::pair<unsigned int, unsigned int>>();
446template <
int dim,
int spacedim>
450 const unsigned int codim)
const
460 if (this->degree < fe_faceq_other->degree)
462 else if (this->degree == fe_faceq_other->degree)
470 if (fe_nothing->is_dominating())
483template <
int dim,
int spacedim>
484std::pair<Table<2, bool>, std::vector<unsigned int>>
488 for (
unsigned int i = 0; i < this->n_dofs_per_cell(); ++i)
489 constant_modes(0, i) =
true;
490 return std::pair<Table<2, bool>, std::vector<unsigned int>>(
491 constant_modes, std::vector<unsigned int>(1, 0));
494template <
int dim,
int spacedim>
498 std::vector<double> & nodal_values)
const
501 this->get_unit_support_points().size());
505 for (
unsigned int i = 0; i < this->n_dofs_per_cell(); ++i)
509 nodal_values[i] = support_point_values[i](0);
515template <
int spacedim>
535template <
int spacedim>
536std::unique_ptr<FiniteElement<1, spacedim>>
539 return std::make_unique<FE_FaceQ<1, spacedim>>(this->degree);
544template <
int spacedim>
551 std::ostringstream namebuf;
553 << this->degree <<
")";
555 return namebuf.str();
560template <
int spacedim>
565 const unsigned int face_no)
const
567 get_subface_interpolation_matrix(source_fe,
569 interpolation_matrix,
575template <
int spacedim>
581 const unsigned int face_no)
const
586 Assert(interpolation_matrix.
n() == this->n_dofs_per_face(face_no),
588 this->n_dofs_per_face(face_no)));
592 interpolation_matrix(0, 0) = 1.;
597template <
int spacedim>
600 const unsigned int face_index)
const
603 return (face_index == shape_index);
608template <
int spacedim>
609std::vector<unsigned int>
612 std::vector<unsigned int> dpo(2, 0U);
619template <
int spacedim>
626template <
int spacedim>
627std::vector<std::pair<unsigned int, unsigned int>>
632 return std::vector<std::pair<unsigned int, unsigned int>>(1,
639template <
int spacedim>
640std::vector<std::pair<unsigned int, unsigned int>>
645 return std::vector<std::pair<unsigned int, unsigned int>>();
650template <
int spacedim>
651std::vector<std::pair<unsigned int, unsigned int>>
654 const unsigned int)
const
657 return std::vector<std::pair<unsigned int, unsigned int>>();
662template <
int spacedim>
663std::pair<Table<2, bool>, std::vector<unsigned int>>
667 for (
unsigned int i = 0; i < this->n_dofs_per_cell(); ++i)
668 constant_modes(0, i) =
true;
669 return std::pair<Table<2, bool>, std::vector<unsigned int>>(
670 constant_modes, std::vector<unsigned int>(1, 0));
675template <
int spacedim>
691template <
int spacedim>
699 const ::internal::FEValuesImplementation::MappingRelatedData<1,
712template <
int spacedim>
716 const unsigned int face,
720 const ::internal::FEValuesImplementation::MappingRelatedData<1,
728 const unsigned int foffset = face;
731 for (
unsigned int k = 0; k < this->n_dofs_per_cell(); ++k)
732 output_data.shape_values(k, 0) = 0.;
733 output_data.shape_values(foffset, 0) = 1;
738template <
int spacedim>
747 const ::internal::FEValuesImplementation::MappingRelatedData<1,
762template <
int dim,
int spacedim>
766 Polynomials::Legendre::generate_complete_basis(degree)),
776template <
int dim,
int spacedim>
777std::unique_ptr<FiniteElement<dim, spacedim>>
780 return std::make_unique<FE_FaceP<dim, spacedim>>(this->degree);
785template <
int dim,
int spacedim>
792 std::ostringstream namebuf;
794 << this->degree <<
")";
796 return namebuf.str();
801template <
int dim,
int spacedim>
804 const unsigned int shape_index,
805 const unsigned int face_index)
const
807 return (face_index == (shape_index / this->n_dofs_per_face(face_index)));
812template <
int dim,
int spacedim>
813std::vector<unsigned int>
816 std::vector<unsigned int> dpo(dim + 1, 0U);
817 dpo[dim - 1] = deg + 1;
818 for (
unsigned int i = 1; i < dim - 1; ++i)
820 dpo[dim - 1] *= deg + 1 + i;
821 dpo[dim - 1] /= i + 1;
828template <
int dim,
int spacedim>
837template <
int dim,
int spacedim>
841 const unsigned int codim)
const
851 if (this->degree < fe_facep_other->degree)
853 else if (this->degree == fe_facep_other->degree)
861 if (fe_nothing->is_dominating())
876template <
int dim,
int spacedim>
881 const unsigned int face_no)
const
883 get_subface_interpolation_matrix(source_fe,
885 interpolation_matrix,
891template <
int dim,
int spacedim>
895 const unsigned int subface,
897 const unsigned int face_no)
const
901 Assert(interpolation_matrix.
n() == this->n_dofs_per_face(face_no),
903 this->n_dofs_per_face(face_no)));
918 this->n_dofs_per_face(face_no) <= source_fe->n_dofs_per_face(face_no),
920 spacedim>::ExcInterpolationNotImplemented()));
925 const QGauss<dim - 1> face_quadrature(source_fe->degree + 1);
930 const double eps = 2e-13 * (this->degree + 1) * (dim - 1);
933 source_fe->n_dofs_per_face(face_no));
935 for (
unsigned int k = 0; k < face_quadrature.size(); ++k)
937 const Point<dim - 1> p =
939 face_quadrature.point(k) :
941 face_quadrature.point(k), subface);
943 for (
unsigned int j = 0; j < source_fe->n_dofs_per_face(face_no); ++j)
944 mass(k, j) = source_fe->poly_space.compute_value(j, p);
954 for (
unsigned int i = 0; i < this->n_dofs_per_face(face_no); ++i)
956 for (
unsigned int k = 0; k < face_quadrature.size(); ++k)
958 const Point<dim - 1> p =
960 face_quadrature.point(k) :
962 face_quadrature.point(k), subface);
963 v_in(k) = this->poly_space.compute_value(i, p);
969 for (
unsigned int j = 0; j < source_fe->n_dofs_per_face(face_no); ++j)
971 double matrix_entry = v_out(j);
976 if (std::fabs(matrix_entry - 1.0) < eps)
978 if (std::fabs(matrix_entry) < eps)
981 interpolation_matrix(j, i) = matrix_entry;
985 else if (
dynamic_cast<const FE_Nothing<dim> *
>(&x_source_fe) !=
nullptr)
993 spacedim>::ExcInterpolationNotImplemented()));
998template <
int dim,
int spacedim>
999std::pair<Table<2, bool>, std::vector<unsigned int>>
1004 constant_modes(0, face * this->n_dofs_per_face(face)) =
true;
1005 return std::pair<Table<2, bool>, std::vector<unsigned int>>(
1006 constant_modes, std::vector<unsigned int>(1, 0));
1011template <
int spacedim>
1018template <
int spacedim>
1025 std::ostringstream namebuf;
1027 << this->degree <<
")";
1029 return namebuf.str();
1035#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
std::vector< std::vector< Point< dim - 1 > > > unit_face_support_points
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 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 UpdateFlags requires_update_flags(const UpdateFlags update_flags) const =0
std::vector< Point< dim > > unit_support_points
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
Abstract base class for mapping classes.
const std::vector< Point< dim > > & get_points() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
@ 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.
static ::ExceptionBase & ExcNotImplemented()
#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)
double least_squares(Vector< number2 > &dst, const Vector< number2 > &src) const
@ 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 Point< dim > child_to_cell_coordinates(const Point< dim > &p, const unsigned int child_index, const RefinementCase< dim > refine_case=RefinementCase< dim >::isotropic_refinement)