22 #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));
50 template <
int dim,
int spacedim>
55 internal::FE_FaceQImplementation::get_QGaussLobatto_points(degree))),
60 std::vector<
bool>(1, true))
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)
120 template <
int dim,
int spacedim>
121 std::unique_ptr<FiniteElement<dim, spacedim>>
124 return std_cxx14::make_unique<FE_FaceQ<dim, spacedim>>(this->degree);
129 template <
int dim,
int spacedim>
136 std::ostringstream namebuf;
138 << this->degree <<
")";
140 return namebuf.str();
145 template <
int dim,
int spacedim>
151 get_subface_interpolation_matrix(source_fe,
153 interpolation_matrix);
158 template <
int dim,
int spacedim>
162 const unsigned int subface,
167 Assert(interpolation_matrix.
n() == this->dofs_per_face,
183 this->dofs_per_face <= source_fe->dofs_per_face,
185 spacedim>::ExcInterpolationNotImplemented()));
189 source_fe->get_unit_face_support_points());
194 const double eps = 2
e-13 * (this->degree + 1) * (dim - 1);
198 for (
unsigned int i = 0; i < source_fe->dofs_per_face; ++i)
200 const Point<dim - 1> p =
202 face_quadrature.point(i) :
204 face_quadrature.point(i), subface);
206 for (
unsigned int j = 0; j < this->dofs_per_face; ++j)
208 double matrix_entry = this->poly_space.compute_value(j, p);
218 interpolation_matrix(i, j) = matrix_entry;
224 for (
unsigned int j = 0; j < source_fe->dofs_per_face; ++j)
228 for (
unsigned int i = 0; i < this->dofs_per_face; ++i)
229 sum += interpolation_matrix(j, i);
234 else if (
dynamic_cast<const FE_Nothing<dim> *
>(&x_source_fe) !=
nullptr)
242 spacedim>::ExcInterpolationNotImplemented()));
247 template <
int dim,
int spacedim>
250 const unsigned int shape_index,
251 const unsigned int face_index)
const
253 return (face_index == (shape_index / this->dofs_per_face));
258 template <
int dim,
int spacedim>
259 std::vector<unsigned int>
262 std::vector<unsigned int> dpo(dim + 1, 0
U);
263 dpo[dim - 1] = deg + 1;
264 for (
unsigned int i = 1; i < dim - 1; ++i)
265 dpo[dim - 1] *= deg + 1;
271 template <
int dim,
int spacedim>
280 template <
int dim,
int spacedim>
281 std::vector<std::pair<unsigned int, unsigned int>>
286 return std::vector<std::pair<unsigned int, unsigned int>>();
291 template <
int dim,
int spacedim>
292 std::vector<std::pair<unsigned int, unsigned int>>
300 return std::vector<std::pair<unsigned int, unsigned int>>();
313 const unsigned int p = this->degree;
314 const unsigned int q = fe_q_other->degree;
316 std::vector<std::pair<unsigned int, unsigned int>> identities;
318 const std::vector<unsigned int> &index_map_inverse =
319 this->poly_space.get_numbering_inverse();
320 const std::vector<unsigned int> &index_map_inverse_other =
321 fe_q_other->poly_space.get_numbering_inverse();
323 for (
unsigned int i = 0; i < p + 1; ++i)
324 for (
unsigned int j = 0; j < q + 1; ++j)
326 this->unit_support_points[index_map_inverse[i]][dim - 1] -
327 fe_q_other->unit_support_points[index_map_inverse_other[j]]
329 identities.emplace_back(i, j);
337 return std::vector<std::pair<unsigned int, unsigned int>>();
348 return std::vector<std::pair<unsigned int, unsigned int>>();
353 return std::vector<std::pair<unsigned int, unsigned int>>();
360 template <
int dim,
int spacedim>
361 std::vector<std::pair<unsigned int, unsigned int>>
369 return std::vector<std::pair<unsigned int, unsigned int>>();
381 const unsigned int p = this->degree;
382 const unsigned int q = fe_q_other->degree;
384 std::vector<std::pair<unsigned int, unsigned int>> identities;
386 const std::vector<unsigned int> &index_map_inverse =
387 this->poly_space.get_numbering_inverse();
388 const std::vector<unsigned int> &index_map_inverse_other =
389 fe_q_other->poly_space.get_numbering_inverse();
391 std::vector<std::pair<unsigned int, unsigned int>> identities_1d;
393 for (
unsigned int i = 0; i < p + 1; ++i)
394 for (
unsigned int j = 0; j < q + 1; ++j)
396 this->unit_support_points[index_map_inverse[i]][dim - 2] -
397 fe_q_other->unit_support_points[index_map_inverse_other[j]]
399 identities_1d.emplace_back(i, j);
401 for (
unsigned int n1 = 0; n1 < identities_1d.size(); ++n1)
402 for (
unsigned int n2 = 0; n2 < identities_1d.size(); ++n2)
403 identities.emplace_back(identities_1d[n1].first * (p + 1) +
404 identities_1d[n2].first,
405 identities_1d[n1].second * (q + 1) +
406 identities_1d[n2].
second);
414 return std::vector<std::pair<unsigned int, unsigned int>>();
425 return std::vector<std::pair<unsigned int, unsigned int>>();
430 return std::vector<std::pair<unsigned int, unsigned int>>();
437 template <
int dim,
int spacedim>
441 const unsigned int codim)
const
451 if (this->degree < fe_faceq_other->degree)
453 else if (this->degree == fe_faceq_other->degree)
461 if (fe_nothing->is_dominating())
474 template <
int dim,
int spacedim>
475 std::pair<Table<2, bool>, std::vector<unsigned int>>
479 for (
unsigned int i = 0; i < this->dofs_per_cell; ++i)
480 constant_modes(0, i) =
true;
481 return std::pair<Table<2, bool>, std::vector<unsigned int>>(
482 constant_modes, std::vector<unsigned int>(1, 0));
485 template <
int dim,
int spacedim>
489 std::vector<double> & nodal_values)
const
492 this->get_unit_support_points().size());
496 for (
unsigned int i = 0; i < this->dofs_per_cell; ++i)
500 nodal_values[i] = support_point_values[i](0);
506 template <
int spacedim>
513 std::vector<
bool>(1, true),
526 template <
int spacedim>
527 std::unique_ptr<FiniteElement<1, spacedim>>
530 return std_cxx14::make_unique<FE_FaceQ<1, spacedim>>(this->
degree);
535 template <
int spacedim>
542 std::ostringstream namebuf;
546 return namebuf.str();
551 template <
int spacedim>
559 interpolation_matrix);
564 template <
int spacedim>
572 Assert(interpolation_matrix.
n() == this->dofs_per_face,
577 interpolation_matrix(0, 0) = 1.;
582 template <
int spacedim>
585 const unsigned int face_index)
const
588 return (face_index == shape_index);
593 template <
int spacedim>
594 std::vector<unsigned int>
597 std::vector<unsigned int> dpo(2, 0
U);
604 template <
int spacedim>
611 template <
int spacedim>
612 std::vector<std::pair<unsigned int, unsigned int>>
617 return std::vector<std::pair<unsigned int, unsigned int>>(1,
624 template <
int spacedim>
625 std::vector<std::pair<unsigned int, unsigned int>>
630 return std::vector<std::pair<unsigned int, unsigned int>>();
635 template <
int spacedim>
636 std::vector<std::pair<unsigned int, unsigned int>>
641 return std::vector<std::pair<unsigned int, unsigned int>>();
646 template <
int spacedim>
647 std::pair<Table<2, bool>, std::vector<unsigned int>>
652 constant_modes(0, i) =
true;
653 return std::pair<Table<2, bool>, std::vector<unsigned int>>(
654 constant_modes, std::vector<unsigned int>(1, 0));
659 template <
int spacedim>
675 template <
int spacedim>
683 const ::internal::FEValuesImplementation::MappingRelatedData<1,
696 template <
int spacedim>
700 const unsigned int face,
704 const ::internal::FEValuesImplementation::MappingRelatedData<1,
712 const unsigned int foffset = face;
716 output_data.shape_values(k, 0) = 0.;
717 output_data.shape_values(foffset, 0) = 1;
722 template <
int spacedim>
731 const ::internal::FEValuesImplementation::MappingRelatedData<1,
746 template <
int dim,
int spacedim>
750 Polynomials::Legendre::generate_complete_basis(degree)),
755 std::vector<
bool>(1, true))
760 template <
int dim,
int spacedim>
761 std::unique_ptr<FiniteElement<dim, spacedim>>
764 return std_cxx14::make_unique<FE_FaceP<dim, spacedim>>(this->degree);
769 template <
int dim,
int spacedim>
776 std::ostringstream namebuf;
778 << this->degree <<
")";
780 return namebuf.str();
785 template <
int dim,
int spacedim>
788 const unsigned int shape_index,
789 const unsigned int face_index)
const
791 return (face_index == (shape_index / this->dofs_per_face));
796 template <
int dim,
int spacedim>
797 std::vector<unsigned int>
800 std::vector<unsigned int> dpo(dim + 1, 0
U);
801 dpo[dim - 1] = deg + 1;
802 for (
unsigned int i = 1; i < dim - 1; ++i)
804 dpo[dim - 1] *= deg + 1 + i;
805 dpo[dim - 1] /= i + 1;
812 template <
int dim,
int spacedim>
821 template <
int dim,
int spacedim>
825 const unsigned int codim)
const
835 if (this->degree < fe_facep_other->degree)
837 else if (this->degree == fe_facep_other->degree)
845 if (fe_nothing->is_dominating())
860 template <
int dim,
int spacedim>
866 get_subface_interpolation_matrix(source_fe,
868 interpolation_matrix);
873 template <
int dim,
int spacedim>
877 const unsigned int subface,
882 Assert(interpolation_matrix.
n() == this->dofs_per_face,
898 this->dofs_per_face <= source_fe->dofs_per_face,
900 spacedim>::ExcInterpolationNotImplemented()));
905 const QGauss<dim - 1> face_quadrature(source_fe->degree + 1);
910 const double eps = 2
e-13 * (this->degree + 1) * (dim - 1);
914 for (
unsigned int k = 0; k < face_quadrature.size(); ++k)
916 const Point<dim - 1> p =
918 face_quadrature.point(k) :
920 face_quadrature.point(k), subface);
922 for (
unsigned int j = 0; j < source_fe->dofs_per_face; ++j)
923 mass(k, j) = source_fe->poly_space.compute_value(j, p);
933 for (
unsigned int i = 0; i < this->dofs_per_face; ++i)
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);
942 v_in(k) = this->poly_space.compute_value(i, p);
948 for (
unsigned int j = 0; j < source_fe->dofs_per_face; ++j)
950 double matrix_entry = v_out(j);
960 interpolation_matrix(j, i) = matrix_entry;
964 else if (
dynamic_cast<const FE_Nothing<dim> *
>(&x_source_fe) !=
nullptr)
972 spacedim>::ExcInterpolationNotImplemented()));
977 template <
int dim,
int spacedim>
978 std::pair<Table<2, bool>, std::vector<unsigned int>>
983 constant_modes(0, face * this->dofs_per_face) =
true;
984 return std::pair<Table<2, bool>, std::vector<unsigned int>>(
985 constant_modes, std::vector<unsigned int>(1, 0));
990 template <
int spacedim>
997 template <
int spacedim>
1004 std::ostringstream namebuf;
1008 return namebuf.str();
1014 #include "fe_face.inst"