45 PolynomialsRaviartThomasNodal(
const unsigned int degree);
66 name()
const override;
74 n_polynomials(
const unsigned int degree);
76 const std::vector<unsigned int> &
77 get_renumbering()
const;
82 virtual std::unique_ptr<TensorPolynomialsBase<dim>>
83 clone()
const override;
92 std::vector<Point<dim>>
93 get_polynomial_support_points()
const;
110 std::vector<unsigned int> lexicographic_to_hierarchic;
116 std::vector<unsigned int> hierarchic_to_lexicographic;
121 std::array<std::vector<unsigned int>, dim> renumber_aniso;
131 std::vector<std::vector<Polynomials::Polynomial<double>>>
132 create_rt_polynomials(
const unsigned int dim,
const unsigned int degree)
134 std::vector<std::vector<Polynomials::Polynomial<double>>> pols(dim);
138 for (
unsigned int d = 1;
d < dim; ++
d)
142 for (
unsigned int d = 1;
d < dim; ++
d)
151 PolynomialsRaviartThomasNodal<dim>::PolynomialsRaviartThomasNodal(
152 const unsigned int degree)
155 , polynomial_space(create_rt_polynomials(dim, degree))
160 const unsigned int n_pols = polynomial_space.n();
161 lexicographic_to_hierarchic =
162 internal::get_lexicographic_numbering_rt_nodal<dim>(degree + 1);
164 hierarchic_to_lexicographic =
170 renumber_aniso[0].resize(n_pols);
171 for (
unsigned int i = 0; i < n_pols; ++i)
172 renumber_aniso[0][i] = i;
176 renumber_aniso[1].resize(n_pols);
177 for (
unsigned int j = 0; j < degree + 2; ++j)
178 for (
unsigned int i = 0; i < degree + 1; ++i)
179 renumber_aniso[1][j * (degree + 1) + i] = j + i * (degree + 2);
184 renumber_aniso[1].resize(n_pols);
185 for (
unsigned int k = 0; k < degree + 1; ++k)
186 for (
unsigned int j = 0; j < degree + 2; ++j)
187 for (
unsigned int i = 0; i < degree + 1; ++i)
188 renumber_aniso[1][(k * (degree + 2) + j) * (degree + 1) + i] =
189 j + k * (degree + 2) + i * (degree + 2) * (degree + 1);
192 renumber_aniso[2].resize(n_pols);
193 for (
unsigned int k = 0; k < degree + 2; ++k)
194 for (
unsigned int j = 0; j < degree + 1; ++j)
195 for (
unsigned int i = 0; i < degree + 1; ++i)
196 renumber_aniso[2][(k * (degree + 1) + j) * (degree + 1) + i] =
197 k + i * (degree + 2) + j * (degree + 2) * (degree + 1);
205 PolynomialsRaviartThomasNodal<dim>::evaluate(
215 Assert(grads.size() == this->n() || grads.size() == 0,
217 Assert(grad_grads.size() == this->n() || grad_grads.size() == 0,
219 Assert(third_derivatives.size() == this->n() ||
220 third_derivatives.size() == 0,
222 Assert(fourth_derivatives.size() == this->n() ||
223 fourth_derivatives.size() == 0,
226 std::vector<double> p_values;
227 std::vector<Tensor<1, dim>> p_grads;
228 std::vector<Tensor<2, dim>> p_grad_grads;
229 std::vector<Tensor<3, dim>> p_third_derivatives;
230 std::vector<Tensor<4, dim>> p_fourth_derivatives;
232 const unsigned int n_sub = polynomial_space.n();
233 p_values.resize((
values.size() == 0) ? 0 : n_sub);
234 p_grads.resize((grads.size() == 0) ? 0 : n_sub);
235 p_grad_grads.resize((grad_grads.size() == 0) ? 0 : n_sub);
236 p_third_derivatives.resize((third_derivatives.size() == 0) ? 0 : n_sub);
237 p_fourth_derivatives.resize((fourth_derivatives.size() == 0) ? 0 : n_sub);
239 for (
unsigned int d = 0;
d < dim; ++
d)
247 for (
unsigned int c = 0; c < dim; ++c)
248 p(c) = unit_point((c + d) % dim);
250 polynomial_space.evaluate(p,
255 p_fourth_derivatives);
257 for (
unsigned int i = 0; i < p_values.size(); ++i)
258 values[lexicographic_to_hierarchic[i + d * n_sub]][d] =
259 p_values[renumber_aniso[d][i]];
261 for (
unsigned int i = 0; i < p_grads.size(); ++i)
262 for (
unsigned int d1 = 0; d1 < dim; ++d1)
263 grads[lexicographic_to_hierarchic[i + d * n_sub]][d]
264 [(d1 + d) % dim] = p_grads[renumber_aniso[
d][i]][d1];
266 for (
unsigned int i = 0; i < p_grad_grads.size(); ++i)
267 for (
unsigned int d1 = 0; d1 < dim; ++d1)
268 for (
unsigned int d2 = 0; d2 < dim; ++d2)
269 grad_grads[lexicographic_to_hierarchic[i + d * n_sub]][d]
270 [(d1 + d) % dim][(d2 +
d) % dim] =
271 p_grad_grads[renumber_aniso[d][i]][d1][d2];
273 for (
unsigned int i = 0; i < p_third_derivatives.size(); ++i)
274 for (
unsigned int d1 = 0; d1 < dim; ++d1)
275 for (
unsigned int d2 = 0; d2 < dim; ++d2)
276 for (
unsigned int d3 = 0; d3 < dim; ++d3)
277 third_derivatives[lexicographic_to_hierarchic[i + d * n_sub]][d]
278 [(d1 + d) % dim][(d2 +
d) % dim]
280 p_third_derivatives[renumber_aniso[
d][i]][d1]
283 for (
unsigned int i = 0; i < p_fourth_derivatives.size(); ++i)
284 for (
unsigned int d1 = 0; d1 < dim; ++d1)
285 for (
unsigned int d2 = 0; d2 < dim; ++d2)
286 for (
unsigned int d3 = 0; d3 < dim; ++d3)
287 for (
unsigned int d4 = 0; d4 < dim; ++d4)
288 fourth_derivatives[lexicographic_to_hierarchic[i + d * n_sub]]
289 [d][(d1 + d) % dim][(d2 +
d) % dim]
290 [(d3 + d) % dim][(d4 +
d) % dim] =
291 p_fourth_derivatives[renumber_aniso[d][i]]
300 PolynomialsRaviartThomasNodal<dim>::name()
const
302 return "PolynomialsRaviartThomasNodal";
309 PolynomialsRaviartThomasNodal<dim>::n_polynomials(
unsigned int degree)
317 std::unique_ptr<TensorPolynomialsBase<dim>>
318 PolynomialsRaviartThomasNodal<dim>::clone()
const
320 return std::make_unique<PolynomialsRaviartThomasNodal<dim>>(*this);
326 std::vector<Point<dim>>
327 PolynomialsRaviartThomasNodal<dim>::get_polynomial_support_points()
const
339 const unsigned int n_sub = polynomial_space.n();
340 std::vector<Point<dim>> points(dim * n_sub);
341 points.resize(n_polynomials(degree));
342 for (
unsigned int d = 0;
d < dim; ++
d)
343 for (
unsigned int i = 0; i < n_sub; ++i)
344 for (
unsigned int c = 0; c < dim; ++c)
345 points[lexicographic_to_hierarchic[i + d * n_sub]][(c + d) % dim] =
346 quad.
point(renumber_aniso[d][i])[c];
358 std::vector<unsigned int>
359 get_rt_dpo_vector(
const unsigned int dim,
const unsigned int degree)
361 std::vector<unsigned int> dpo(dim + 1);
364 unsigned int dofs_per_face = 1;
365 for (
unsigned int d = 1;
d < dim; ++
d)
366 dofs_per_face *= (degree + 1);
368 dpo[dim - 1] = dofs_per_face;
369 dpo[dim] = dim * degree * dofs_per_face;
378 std::vector<unsigned int>
382 std::vector<unsigned int> lexicographic_numbering;
384 for (
unsigned int j = 0; j < n_dofs_face; j++)
386 lexicographic_numbering.push_back(j);
387 for (
unsigned int i = n_dofs_face * 2 * dim;
388 i < n_dofs_face * 2 * dim + degree - 1;
390 lexicographic_numbering.push_back(i + j * (degree - 1));
391 lexicographic_numbering.push_back(n_dofs_face + j);
395 unsigned int layers = (dim == 3) ? degree : 1;
396 for (
unsigned int k = 0; k < layers; k++)
398 unsigned int k_add = k * degree;
399 for (
unsigned int j = n_dofs_face * 2; j < n_dofs_face * 2 + degree;
401 lexicographic_numbering.push_back(j + k_add);
403 for (
unsigned int i = n_dofs_face * (2 * dim + (degree - 1));
404 i < n_dofs_face * (2 * dim + (degree - 1)) + (degree - 1) * degree;
407 lexicographic_numbering.push_back(i + k_add * (degree - 1));
409 for (
unsigned int j = n_dofs_face * 3; j < n_dofs_face * 3 + degree;
411 lexicographic_numbering.push_back(j + k_add);
417 for (
unsigned int i = 4 * n_dofs_face; i < 5 * n_dofs_face; i++)
418 lexicographic_numbering.push_back(i);
419 for (
unsigned int i = 6 * n_dofs_face + n_dofs_face * 2 * (degree - 1);
420 i < 6 * n_dofs_face + n_dofs_face * 3 * (degree - 1);
422 lexicographic_numbering.push_back(i);
423 for (
unsigned int i = 5 * n_dofs_face; i < 6 * n_dofs_face; i++)
424 lexicographic_numbering.push_back(i);
427 return lexicographic_numbering;
437 :
FE_PolyTensor<dim>(PolynomialsRaviartThomasNodal<dim>(degree),
444 PolynomialsRaviartThomasNodal<dim>::n_polynomials(
455 PolynomialsRaviartThomasNodal<dim>(
degree).get_polynomial_support_points();
457 this->n_dofs_per_cell());
459 const unsigned int face_no = 0;
466 for (
unsigned int i = 0; i < GeometryInfo<dim>::max_children_per_face; ++i)
469 FETools::compute_face_embedding_matrices<dim, double>(*
this,
476 unsigned int target_row = 0;
477 for (
unsigned int d = 0; d < GeometryInfo<dim>::max_children_per_face; ++d)
478 for (
unsigned int i = 0; i < face_embeddings[d].
m(); ++i)
480 for (
unsigned int j = 0; j < face_embeddings[d].
n(); ++j)
502 return "FE_RaviartThomasNodal<" + std::to_string(dim) +
">(" +
503 std::to_string(this->degree - 1) +
")";
508std::unique_ptr<FiniteElement<dim, dim>>
511 return std::make_unique<FE_RaviartThomasNodal<dim>>(*this);
524 dim>::initialize_quad_dof_index_permutation_and_sign_change()
530 const unsigned int n = this->degree;
531 const unsigned int face_no = 0;
533 for (
unsigned int local = 0; local < this->n_dofs_per_quad(face_no); ++local)
537 unsigned int i = local % n, j = local / n;
540 this->adjust_quad_dof_index_for_face_orientation_table[face_no](local,
544 this->adjust_quad_dof_index_for_face_orientation_table[face_no](local,
546 i + (n - 1 - j) * n - local;
548 this->adjust_quad_dof_index_for_face_orientation_table[face_no](local,
550 (n - 1 - j) + (n - 1 - i) * n - local;
552 this->adjust_quad_dof_index_for_face_orientation_table[face_no](local,
554 (n - 1 - i) + j * n - local;
556 this->adjust_quad_dof_index_for_face_orientation_table[face_no](local,
559 this->adjust_quad_dof_index_for_face_orientation_table[face_no](local,
561 j + (n - 1 - i) * n - local;
563 this->adjust_quad_dof_index_for_face_orientation_table[face_no](local,
565 (n - 1 - i) + (n - 1 - j) * n - local;
567 this->adjust_quad_dof_index_for_face_orientation_table[face_no](local,
569 (n - 1 - j) + i * n - local;
572 for (
unsigned int i = 0; i < 4; ++i)
573 this->adjust_quad_dof_sign_for_face_orientation_table[face_no](local,
583 const unsigned int shape_index,
584 const unsigned int face_index)
const
591 const unsigned int support_face = shape_index / this->n_dofs_per_face();
609 std::vector<double> & nodal_values)
const
611 Assert(support_point_values.size() == this->generalized_support_points.size(),
613 this->generalized_support_points.size()));
614 Assert(nodal_values.size() == this->n_dofs_per_cell(),
616 Assert(support_point_values[0].size() == this->n_components(),
618 this->n_components()));
622 unsigned int fbase = 0;
624 for (; f < GeometryInfo<dim>::faces_per_cell;
625 ++f, fbase += this->n_dofs_per_face(f))
627 for (
unsigned int i = 0; i < this->n_dofs_per_face(f); ++i)
629 nodal_values[fbase + i] = support_point_values[fbase + i](
635 const unsigned int istep = (this->n_dofs_per_cell() - fbase) / dim;
639 while (fbase < this->n_dofs_per_cell())
641 for (
unsigned int i = 0; i < istep; ++i)
643 nodal_values[fbase + i] = support_point_values[fbase + i](f);
668std::vector<std::pair<unsigned int, unsigned int>>
676 return std::vector<std::pair<unsigned int, unsigned int>>();
678 return std::vector<std::pair<unsigned int, unsigned int>>();
682 return std::vector<std::pair<unsigned int, unsigned int>>();
689std::vector<std::pair<unsigned int, unsigned int>>
700 return std::vector<std::pair<unsigned int, unsigned int>>();
715 const unsigned int p = this->degree - 1;
716 const unsigned int q = fe_q_other->degree - 1;
718 std::vector<std::pair<unsigned int, unsigned int>> identities;
721 for (
unsigned int i = 0; i < p + 1; ++i)
722 identities.emplace_back(i, i);
724 else if (p % 2 == 0 && q % 2 == 0)
725 identities.emplace_back(p / 2, q / 2);
733 return std::vector<std::pair<unsigned int, unsigned int>>();
738 return std::vector<std::pair<unsigned int, unsigned int>>();
744std::vector<std::pair<unsigned int, unsigned int>>
747 const unsigned int face_no)
const
756 return std::vector<std::pair<unsigned int, unsigned int>>();
759 const unsigned int p = this->n_dofs_per_quad(face_no);
762 const unsigned int q = fe_q_other->n_dofs_per_quad(0);
764 std::vector<std::pair<unsigned int, unsigned int>> identities;
767 for (
unsigned int i = 0; i < p; ++i)
768 identities.emplace_back(i, i);
770 else if (p % 2 != 0 && q % 2 != 0)
771 identities.emplace_back(p / 2, q / 2);
779 return std::vector<std::pair<unsigned int, unsigned int>>();
784 return std::vector<std::pair<unsigned int, unsigned int>>();
793 const unsigned int codim)
const
803 if (this->degree < fe_rt_nodal_other->degree)
805 else if (this->degree == fe_rt_nodal_other->degree)
813 if (fe_nothing->is_dominating())
833 const unsigned int)
const
845 const unsigned int)
const
857 const unsigned int face_no)
const
864 &x_source_fe) !=
nullptr),
867 Assert(interpolation_matrix.
n() == this->n_dofs_per_face(face_no),
869 this->n_dofs_per_face(face_no)));
895 double eps = 2e-13 * this->degree * (dim - 1);
906 const Point<dim> &p = face_projection.point(i);
908 for (
unsigned int j = 0; j < this->n_dofs_per_face(face_no); ++j)
910 double matrix_entry =
911 this->shape_value_component(this->face_to_cell_index(j, 0), p, 0);
916 if (std::fabs(matrix_entry - 1.0) < eps)
918 if (std::fabs(matrix_entry) < eps)
921 interpolation_matrix(i, j) = matrix_entry;
932 for (
unsigned int i = 0; i < this->n_dofs_per_face(face_no); ++i)
933 sum += interpolation_matrix(j, i);
935 Assert(std::fabs(sum - 1) < 2e-13 * this->degree * (dim - 1),
946 const unsigned int subface,
948 const unsigned int face_no)
const
954 &x_source_fe) !=
nullptr),
957 Assert(interpolation_matrix.
n() == this->n_dofs_per_face(face_no),
959 this->n_dofs_per_face(face_no)));
985 double eps = 2e-13 * this->degree * (dim - 1);
997 const Point<dim> &p = subface_projection.point(i);
999 for (
unsigned int j = 0; j < this->n_dofs_per_face(face_no); ++j)
1001 double matrix_entry =
1002 this->shape_value_component(this->face_to_cell_index(j, 0), p, 0);
1007 if (std::fabs(matrix_entry - 1.0) < eps)
1009 if (std::fabs(matrix_entry) < eps)
1012 interpolation_matrix(i, j) = matrix_entry;
1023 for (
unsigned int i = 0; i < this->n_dofs_per_face(face_no); ++i)
1024 sum += interpolation_matrix(j, i);
1026 Assert(std::fabs(sum - 1) < 2e-13 * this->degree * (dim - 1),
1037 const unsigned int child,
1044 "Prolongation matrices are only available for refined cells!"));
1048 if (this->prolongation[refinement_case - 1][child].n() == 0)
1050 std::lock_guard<std::mutex> lock(this->mutex);
1053 if (this->prolongation[refinement_case - 1][child].n() ==
1054 this->n_dofs_per_cell())
1055 return this->prolongation[refinement_case - 1][child];
1063 std::vector<std::vector<FullMatrix<double>>> isotropic_matrices(
1065 isotropic_matrices.back().resize(
1068 this->n_dofs_per_cell()));
1071 isotropic_matrices.back());
1087 return this->prolongation[refinement_case - 1][child];
1095 const unsigned int child,
1102 "Restriction matrices are only available for refined cells!"));
1106 if (this->restriction[refinement_case - 1][child].n() == 0)
1108 std::lock_guard<std::mutex> lock(this->mutex);
1111 if (this->restriction[refinement_case - 1][child].n() ==
1112 this->n_dofs_per_cell())
1113 return this->restriction[refinement_case - 1][child];
1121 std::vector<std::vector<FullMatrix<double>>> isotropic_matrices(
1123 isotropic_matrices.back().resize(
1126 this->n_dofs_per_cell()));
1128 this_nonconst.
restriction[refinement_case - 1].swap(
1129 isotropic_matrices.back());
1145 return this->restriction[refinement_case - 1][child];
1151#include "fe_raviart_thomas_nodal.inst"
std::vector< MappingKind > mapping_kind
virtual std::unique_ptr< FiniteElement< dim, dim > > clone() const override
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 bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const override
virtual const FullMatrix< double > & get_restriction_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const override
void initialize_quad_dof_index_permutation_and_sign_change()
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities(const FiniteElement< dim > &fe_other) const override
virtual void get_face_interpolation_matrix(const FiniteElement< dim > &source, FullMatrix< double > &matrix, const unsigned int face_no=0) const override
virtual std::string get_name() const override
virtual bool hp_constraints_are_implemented() const override
virtual void get_subface_interpolation_matrix(const FiniteElement< dim > &source, const unsigned int subface, FullMatrix< double > &matrix, const unsigned int face_no=0) const override
virtual FiniteElementDomination::Domination compare_for_domination(const FiniteElement< dim > &fe_other, const unsigned int codim=0) const override final
virtual const FullMatrix< double > & get_prolongation_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const override
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities(const FiniteElement< dim > &fe_other) const override
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_quad_dof_identities(const FiniteElement< dim > &fe_other, const unsigned int face_no=0) const override
FE_RaviartThomasNodal(const unsigned int p)
const unsigned int degree
unsigned int n_dofs_per_face(unsigned int face_no=0, unsigned int child=0) const
virtual std::string get_name() const =0
std::vector< std::vector< FullMatrix< double > > > restriction
void reinit_restriction_and_prolongation_matrices(const bool isotropic_restriction_only=false, const bool isotropic_prolongation_only=false)
std::vector< std::vector< Point< dim - 1 > > > generalized_face_support_points
FullMatrix< double > interface_constraints
std::vector< Point< dim > > generalized_support_points
std::vector< std::vector< FullMatrix< double > > > prolongation
static void project_to_subface(const SubQuadrature &quadrature, const unsigned int face_no, const unsigned int subface_no, std::vector< Point< dim > > &q_points, const RefinementCase< dim - 1 > &ref_case=RefinementCase< dim - 1 >::isotropic_refinement)
static void project_to_face(const SubQuadrature &quadrature, const unsigned int face_no, std::vector< Point< dim > > &q_points)
const Point< dim > & point(const unsigned int i) const
virtual std::string name() const =0
virtual void evaluate(const Point< dim > &unit_point, std::vector< Tensor< 1, dim > > &values, std::vector< Tensor< 2, dim > > &grads, std::vector< Tensor< 3, dim > > &grad_grads, std::vector< Tensor< 4, dim > > &third_derivatives, std::vector< Tensor< 5, dim > > &fourth_derivatives) const =0
virtual std::unique_ptr< TensorPolynomialsBase< dim > > clone() const =0
unsigned int degree() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
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)
@ either_element_can_dominate
@ other_element_dominates
@ neither_element_dominates
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)
constexpr T pow(const T base, const int iexp)
std::vector< Integer > invert_permutation(const std::vector< Integer > &permutation)
std::vector< unsigned int > get_lexicographic_numbering_rt_nodal(const unsigned int degree)