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 = 2
e-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);
222 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);
238 else if (
dynamic_cast<const FE_Nothing<dim> *
>(&x_source_fe) !=
nullptr)
246 spacedim>::ExcInterpolationNotImplemented()));
251template <
int dim,
int spacedim>
254 const unsigned int shape_index,
255 const unsigned int face_index)
const
257 return (face_index == (shape_index / this->n_dofs_per_face(face_index)));
262template <
int dim,
int spacedim>
263std::vector<unsigned int>
266 std::vector<unsigned int> dpo(dim + 1, 0
U);
267 dpo[dim - 1] = deg + 1;
268 for (
unsigned int i = 1; i < dim - 1; ++i)
269 dpo[dim - 1] *= deg + 1;
275template <
int dim,
int spacedim>
284template <
int dim,
int spacedim>
285std::vector<std::pair<unsigned int, unsigned int>>
290 return std::vector<std::pair<unsigned int, unsigned int>>();
295template <
int dim,
int spacedim>
296std::vector<std::pair<unsigned int, unsigned int>>
304 return std::vector<std::pair<unsigned int, unsigned int>>();
317 const unsigned int p = this->degree;
318 const unsigned int q = fe_q_other->degree;
320 std::vector<std::pair<unsigned int, unsigned int>> identities;
322 const std::vector<unsigned int> &index_map_inverse =
323 this->poly_space.get_numbering_inverse();
324 const std::vector<unsigned int> &index_map_inverse_other =
325 fe_q_other->poly_space.get_numbering_inverse();
327 for (
unsigned int i = 0; i < p + 1; ++i)
328 for (
unsigned int j = 0; j < q + 1; ++j)
331 fe_q_other->unit_support_points[index_map_inverse_other[j]]
333 identities.emplace_back(i, j);
341 return std::vector<std::pair<unsigned int, unsigned int>>();
353 return std::vector<std::pair<unsigned int, unsigned int>>();
358 return std::vector<std::pair<unsigned int, unsigned int>>();
365template <
int dim,
int spacedim>
366std::vector<std::pair<unsigned int, unsigned int>>
369 const unsigned int)
const
375 return std::vector<std::pair<unsigned int, unsigned int>>();
387 const unsigned int p = this->degree;
388 const unsigned int q = fe_q_other->degree;
390 std::vector<std::pair<unsigned int, unsigned int>> identities;
392 const std::vector<unsigned int> &index_map_inverse =
393 this->poly_space.get_numbering_inverse();
394 const std::vector<unsigned int> &index_map_inverse_other =
395 fe_q_other->poly_space.get_numbering_inverse();
397 std::vector<std::pair<unsigned int, unsigned int>> identities_1d;
399 for (
unsigned int i = 0; i < p + 1; ++i)
400 for (
unsigned int j = 0; j < q + 1; ++j)
403 fe_q_other->unit_support_points[index_map_inverse_other[j]]
405 identities_1d.emplace_back(i, j);
407 for (
unsigned int n1 = 0; n1 < identities_1d.size(); ++n1)
408 for (
unsigned int n2 = 0; n2 < identities_1d.size(); ++n2)
409 identities.emplace_back(identities_1d[n1].first * (p + 1) +
410 identities_1d[n2].first,
411 identities_1d[n1].second * (q + 1) +
412 identities_1d[n2].
second);
420 return std::vector<std::pair<unsigned int, unsigned int>>();
432 return std::vector<std::pair<unsigned int, unsigned int>>();
437 return std::vector<std::pair<unsigned int, unsigned int>>();
444template <
int dim,
int spacedim>
448 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;
553 return namebuf.str();
558template <
int spacedim>
563 const unsigned int face_no)
const
567 interpolation_matrix,
573template <
int spacedim>
579 const unsigned int face_no)
const
584 Assert(interpolation_matrix.
n() == this->n_dofs_per_face(face_no),
586 this->n_dofs_per_face(face_no)));
590 interpolation_matrix(0, 0) = 1.;
595template <
int spacedim>
598 const unsigned int face_index)
const
601 return (face_index == shape_index);
606template <
int spacedim>
607std::vector<unsigned int>
610 std::vector<unsigned int> dpo(2, 0
U);
617template <
int spacedim>
624template <
int spacedim>
625std::vector<std::pair<unsigned int, unsigned int>>
630 return std::vector<std::pair<unsigned int, unsigned int>>(1,
637template <
int spacedim>
638std::vector<std::pair<unsigned int, unsigned int>>
643 return std::vector<std::pair<unsigned int, unsigned int>>();
648template <
int spacedim>
649std::vector<std::pair<unsigned int, unsigned int>>
652 const unsigned int)
const
655 return std::vector<std::pair<unsigned int, unsigned int>>();
660template <
int spacedim>
661std::pair<Table<2, bool>, std::vector<unsigned int>>
666 constant_modes(0, i) =
true;
667 return std::pair<Table<2, bool>, std::vector<unsigned int>>(
668 constant_modes, std::vector<unsigned int>(1, 0));
673template <
int spacedim>
689template <
int spacedim>
697 const ::internal::FEValuesImplementation::MappingRelatedData<1,
710template <
int spacedim>
714 const unsigned int face,
718 const ::internal::FEValuesImplementation::MappingRelatedData<1,
726 const unsigned int foffset = face;
730 output_data.shape_values(k, 0) = 0.;
731 output_data.shape_values(foffset, 0) = 1;
736template <
int spacedim>
745 const ::internal::FEValuesImplementation::MappingRelatedData<1,
760template <
int dim,
int spacedim>
764 Polynomials::Legendre::generate_complete_basis(degree)),
774template <
int dim,
int spacedim>
775std::unique_ptr<FiniteElement<dim, spacedim>>
778 return std::make_unique<FE_FaceP<dim, spacedim>>(this->degree);
783template <
int dim,
int spacedim>
790 std::ostringstream namebuf;
792 << this->degree <<
")";
794 return namebuf.str();
799template <
int dim,
int spacedim>
802 const unsigned int shape_index,
803 const unsigned int face_index)
const
805 return (face_index == (shape_index / this->n_dofs_per_face(face_index)));
810template <
int dim,
int spacedim>
811std::vector<unsigned int>
814 std::vector<unsigned int> dpo(dim + 1, 0
U);
815 dpo[dim - 1] = deg + 1;
816 for (
unsigned int i = 1; i < dim - 1; ++i)
818 dpo[dim - 1] *= deg + 1 + i;
819 dpo[dim - 1] /= i + 1;
826template <
int dim,
int spacedim>
835template <
int dim,
int spacedim>
839 const unsigned int codim)
const
849 if (this->degree < fe_facep_other->degree)
851 else if (this->degree == fe_facep_other->degree)
859 if (fe_nothing->is_dominating())
874template <
int dim,
int spacedim>
879 const unsigned int face_no)
const
881 get_subface_interpolation_matrix(source_fe,
883 interpolation_matrix,
889template <
int dim,
int spacedim>
893 const unsigned int subface,
895 const unsigned int face_no)
const
899 Assert(interpolation_matrix.
n() == this->n_dofs_per_face(face_no),
901 this->n_dofs_per_face(face_no)));
916 this->n_dofs_per_face(face_no) <= source_fe->n_dofs_per_face(face_no),
918 spacedim>::ExcInterpolationNotImplemented()));
923 const QGauss<dim - 1> face_quadrature(source_fe->degree + 1);
928 const double eps = 2
e-13 * (this->degree + 1) * (dim - 1);
931 source_fe->n_dofs_per_face(face_no));
933 for (
unsigned int k = 0; k < face_quadrature.size(); ++k)
935 const Point<dim - 1> p =
937 face_quadrature.point(k) :
939 face_quadrature.point(k), subface);
941 for (
unsigned int j = 0; j < source_fe->n_dofs_per_face(face_no); ++j)
942 mass(k, j) = source_fe->poly_space.compute_value(j, p);
952 for (
unsigned int i = 0; i < this->n_dofs_per_face(face_no); ++i)
954 for (
unsigned int k = 0; k < face_quadrature.size(); ++k)
956 const Point<dim - 1> p =
958 face_quadrature.point(k) :
960 face_quadrature.point(k), subface);
961 v_in(k) = this->poly_space.compute_value(i, p);
967 for (
unsigned int j = 0; j < source_fe->n_dofs_per_face(face_no); ++j)
969 double matrix_entry = v_out(j);
979 interpolation_matrix(j, i) = matrix_entry;
983 else if (
dynamic_cast<const FE_Nothing<dim> *
>(&x_source_fe) !=
nullptr)
991 spacedim>::ExcInterpolationNotImplemented()));
996template <
int dim,
int spacedim>
997std::pair<Table<2, bool>, std::vector<unsigned int>>
1002 constant_modes(0, face * this->n_dofs_per_face(face)) =
true;
1003 return std::pair<Table<2, bool>, std::vector<unsigned int>>(
1004 constant_modes, std::vector<unsigned int>(1, 0));
1009template <
int spacedim>
1016template <
int spacedim>
1023 std::ostringstream namebuf;
1027 return namebuf.str();
1033#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
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
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
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
const unsigned int degree
unsigned int n_dofs_per_cell() const
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
std::vector< Point< dim > > unit_support_points
Abstract base class for mapping classes.
const std::vector< Point< dim > > & get_points() const
SymmetricTensor< rank, dim, Number > sum(const SymmetricTensor< rank, dim, Number > &local, const MPI_Comm &mpi_communicator)
#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
Expression fabs(const Expression &x)
std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
@ either_element_can_dominate
@ other_element_dominates
@ neither_element_dominates
void L2(Vector< number > &result, const FEValuesBase< dim > &fe, const std::vector< double > &input, const double factor=1.)
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
std::vector< Polynomial< double > > generate_complete_Lagrange_basis(const std::vector< Point< 1 > > &points)
std::string dim_string(const int dim, const int spacedim)
std::vector< Point< dim > > unit_support_points(const std::vector< Point< 1 > > &line_support_points, const std::vector< unsigned int > &renumbering)
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)