92 const unsigned int face_no = 0;
97 for (
unsigned int i = 0; i < GeometryInfo<dim>::max_children_per_face; ++i)
100 FETools::compute_face_embedding_matrices<dim, double>(*
this,
107 unsigned int target_row = 0;
108 for (
unsigned int d = 0; d < GeometryInfo<dim>::max_children_per_face; ++d)
109 for (
unsigned int i = 0; i < face_embeddings[d].
m(); ++i)
111 for (
unsigned int j = 0; j < face_embeddings[d].
n(); ++j)
138 std::ostringstream namebuf;
139 namebuf <<
"FE_RaviartThomas<" << dim <<
">(" << this->degree - 1 <<
")";
141 return namebuf.str();
146std::unique_ptr<FiniteElement<dim, dim>>
149 return std::make_unique<FE_RaviartThomas<dim>>(*this);
163 const unsigned int n_interior_points = (deg > 0) ? cell_quadrature.
size() : 0;
168 const unsigned int face_no = 0;
170 unsigned int n_face_points = (dim > 1) ? 1 : 0;
172 for (
unsigned int d = 1; d < dim; ++d)
173 n_face_points *= deg + 1;
176 this->generalized_support_points.resize(
178 this->generalized_face_support_points[face_no].resize(n_face_points);
181 unsigned int current = 0;
185 QGauss<dim - 1> face_points(deg + 1);
189 boundary_weights.reinit(n_face_points, legendre.n());
191 for (
unsigned int k = 0; k < n_face_points; ++k)
193 this->generalized_face_support_points[face_no][k] =
194 face_points.point(k);
198 for (
unsigned int i = 0; i < legendre.n(); ++i)
200 boundary_weights(k, i) =
201 face_points.weight(k) *
202 legendre.compute_value(i, face_points.point(k));
209 for (; current < GeometryInfo<dim>::faces_per_cell * n_face_points;
214 this->generalized_support_points[current] = faces.
point(
217 this->reference_cell(), 0,
true,
false,
false, n_face_points));
225 std::unique_ptr<AnisotropicPolynomials<dim>> polynomials[dim];
226 for (
unsigned int dd = 0; dd < dim; ++dd)
228 std::vector<std::vector<Polynomials::Polynomial<double>>> poly(dim);
229 for (
unsigned int d = 0; d < dim; ++d)
233 polynomials[dd] = std::make_unique<AnisotropicPolynomials<dim>>(poly);
236 interior_weights.reinit(
239 for (
unsigned int k = 0; k < cell_quadrature.
size(); ++k)
241 this->generalized_support_points[current++] = cell_quadrature.
point(k);
242 for (
unsigned int i = 0; i < polynomials[0]->n(); ++i)
243 for (
unsigned int d = 0; d < dim; ++d)
244 interior_weights(k, i, d) =
245 cell_quadrature.
weight(k) *
246 polynomials[d]->compute_value(i, cell_quadrature.
point(k));
249 Assert(current == this->generalized_support_points.size(),
270 const unsigned int face_no = 0;
273 this->adjust_quad_dof_index_for_face_orientation_table[0].n_elements() ==
274 this->reference_cell().n_face_orientations(face_no) *
275 this->n_dofs_per_quad(face_no),
279 const unsigned int n = this->tensor_degree();
333 for (
unsigned int local_face_dof = 0;
334 local_face_dof < this->n_dofs_per_quad(face_no);
338 unsigned int i = local_face_dof % n;
339 unsigned int j = local_face_dof / n;
344 for (
unsigned int case_no = 0; case_no < 8; ++case_no)
351 if (((!face_orientation) && (!face_rotation)) ||
352 ((face_orientation) && (face_rotation)))
356 this->adjust_quad_dof_index_for_face_orientation_table[face_no](
357 local_face_dof, case_no) = j + i * n - local_face_dof;
362 this->adjust_quad_dof_index_for_face_orientation_table[face_no](
363 local_face_dof, case_no) = 0;
367 const unsigned int new_local_face_dof =
369 this->adjust_quad_dof_index_for_face_orientation_table[face_no](
370 local_face_dof, case_no);
372 i = new_local_face_dof % n;
373 j = new_local_face_dof / n;
380 if (!face_flip && face_rotation)
382 this->adjust_quad_dof_sign_for_face_orientation_table[face_no](
383 local_face_dof, case_no) = ((j % 2) == 1);
386 else if (face_flip && !face_rotation)
390 this->adjust_quad_dof_sign_for_face_orientation_table[face_no](
391 local_face_dof, case_no) = ((j % 2) == 1) != ((i % 2) == 1);
394 else if (face_flip && face_rotation)
396 this->adjust_quad_dof_sign_for_face_orientation_table[face_no](
397 local_face_dof, case_no) = ((i % 2) == 1);
406 if (!face_orientation)
407 this->adjust_quad_dof_sign_for_face_orientation_table[face_no](
408 local_face_dof, case_no) =
409 !this->adjust_quad_dof_sign_for_face_orientation_table[face_no](
410 local_face_dof, case_no);
423 for (
unsigned int i = 0; i < GeometryInfo<1>::max_children_per_cell; ++i)
424 this->restriction[0][i].reinit(0, 0);
445 QGauss<dim - 1> q_base(this->degree);
446 const unsigned int n_face_points = q_base.size();
462 for (
unsigned int k = 0; k < q_face.
size(); ++k)
463 for (
unsigned int i = 0; i < this->n_dofs_per_cell(); ++i)
464 cached_values_on_face(i, k) = this->shape_value_component(
467 for (
unsigned int sub = 0; sub < GeometryInfo<dim>::max_children_per_face;
475 this->reference_cell(), q_base, face, sub);
491 for (
unsigned int k = 0; k < n_face_points; ++k)
492 for (
unsigned int i_child = 0; i_child < this->n_dofs_per_cell();
494 for (
unsigned int i_face = 0;
495 i_face < this->n_dofs_per_face(face);
503 this->restriction[iso][child](
504 face * this->n_dofs_per_face(face) + i_face, i_child) +=
505 Utilities::fixed_power<dim - 1>(.5) * q_sub.
weight(k) *
506 cached_values_on_face(i_child, k) *
507 this->shape_value_component(
508 face * this->n_dofs_per_face(face) + i_face,
515 if (this->degree == 1)
520 std::unique_ptr<AnisotropicPolynomials<dim>> polynomials[dim];
521 for (
unsigned int dd = 0; dd < dim; ++dd)
523 std::vector<std::vector<Polynomials::Polynomial<double>>> poly(dim);
524 for (
unsigned int d = 0; d < dim; ++d)
530 polynomials[dd] = std::make_unique<AnisotropicPolynomials<dim>>(poly);
536 const unsigned int face_no = 0;
539 const unsigned int start_cell_dofs =
548 for (
unsigned int k = 0; k < q_cell.
size(); ++k)
549 for (
unsigned int i = 0; i < this->n_dofs_per_cell(); ++i)
550 for (
unsigned int d = 0; d < dim; ++d)
551 cached_values_on_cell(i, k, d) =
552 this->shape_value_component(i, q_cell.
point(k), d);
554 for (
unsigned int child = 0; child < GeometryInfo<dim>::max_children_per_cell;
562 for (
unsigned int k = 0; k < q_sub.
size(); ++k)
563 for (
unsigned int i_child = 0; i_child < this->n_dofs_per_cell();
565 for (
unsigned int d = 0; d < dim; ++d)
566 for (
unsigned int i_weight = 0; i_weight < polynomials[d]->n();
569 this->restriction[iso][child](start_cell_dofs + i_weight * dim +
572 q_sub.
weight(k) * cached_values_on_cell(i_child, k, d) *
573 polynomials[d]->compute_value(i_weight, q_sub.
point(k));
581std::vector<unsigned int>
586 unsigned int dofs_per_face = 1;
587 for (
unsigned int d = 1; d < dim; ++d)
588 dofs_per_face *= deg + 1;
591 const unsigned int interior_dofs = dim * deg * dofs_per_face;
593 std::vector<unsigned int> dpo(dim + 1);
594 dpo[dim - 1] = dofs_per_face;
595 dpo[dim] = interior_dofs;
603std::pair<Table<2, bool>, std::vector<unsigned int>>
607 for (
unsigned int d = 0; d < dim; ++d)
608 for (
unsigned int i = 0; i < this->n_dofs_per_cell(); ++i)
609 constant_modes(d, i) =
true;
610 std::vector<unsigned int> components;
611 for (
unsigned int d = 0; d < dim; ++d)
612 components.push_back(d);
613 return std::pair<Table<2, bool>, std::vector<unsigned int>>(constant_modes,
627 const unsigned int face_index)
const
635 switch (this->degree)
648 return (face_index !=
670 std::vector<double> & nodal_values)
const
672 Assert(support_point_values.size() == this->generalized_support_points.size(),
674 this->generalized_support_points.size()));
675 Assert(nodal_values.size() == this->n_dofs_per_cell(),
677 Assert(support_point_values[0].size() == this->n_components(),
679 this->n_components()));
681 std::fill(nodal_values.begin(), nodal_values.end(), 0.);
683 const unsigned int n_face_points = boundary_weights.size(0);
685 for (
unsigned int k = 0; k < n_face_points; ++k)
686 for (
unsigned int i = 0; i < boundary_weights.size(1); ++i)
688 nodal_values[i + face * this->n_dofs_per_face(face)] +=
689 boundary_weights(k, i) *
690 support_point_values[face * n_face_points + k](
697 const unsigned int face_no = 0;
699 const unsigned int start_cell_dofs =
701 const unsigned int start_cell_points =
704 for (
unsigned int k = 0; k < interior_weights.size(0); ++k)
705 for (
unsigned int i = 0; i < interior_weights.size(1); ++i)
706 for (
unsigned int d = 0; d < dim; ++d)
707 nodal_values[start_cell_dofs + i * dim + d] +=
708 interior_weights(k, i, d) *
709 support_point_values[k + start_cell_points](d);
725#include "fe_raviart_thomas.inst"
FullMatrix< double > inverse_node_matrix
std::vector< MappingKind > mapping_kind
virtual std::size_t memory_consumption() 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
void initialize_restriction()
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const override
void initialize_quad_dof_index_permutation_and_sign_change()
virtual std::unique_ptr< FiniteElement< dim, dim > > clone() const override
friend class FE_RaviartThomas
void initialize_support_points(const unsigned int rt_degree)
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
virtual std::string get_name() const override
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
void reinit_restriction_and_prolongation_matrices(const bool isotropic_restriction_only=false, const bool isotropic_prolongation_only=false)
FullMatrix< double > interface_constraints
std::vector< std::vector< FullMatrix< double > > > prolongation
void invert(const FullMatrix< number2 > &M)
static std::vector< Polynomial< double > > generate_complete_basis(const unsigned int degree)
static Quadrature< dim > project_to_child(const ReferenceCell &reference_cell, const Quadrature< dim > &quadrature, const unsigned int child_no)
static Quadrature< dim > project_to_all_faces(const ReferenceCell &reference_cell, const hp::QCollection< dim - 1 > &quadrature)
static void project_to_subface(const ReferenceCell &reference_cell, 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 ReferenceCell &reference_cell, const SubQuadrature &quadrature, const unsigned int face_no, std::vector< Point< dim > > &q_points)
const Point< dim > & point(const unsigned int i) const
double weight(const unsigned int i) const
unsigned int size() 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)
bool get_bit(const unsigned char number, const unsigned int n)
static unsigned int child_cell_on_face(const RefinementCase< dim > &ref_case, const unsigned int face, const unsigned int subface, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false, const RefinementCase< dim - 1 > &face_refinement_case=RefinementCase< dim - 1 >::isotropic_refinement)
static std_cxx20::ranges::iota_view< unsigned int, unsigned int > face_indices()