89 const unsigned int face_no = 0;
94 for (
unsigned int i = 0; i < GeometryInfo<dim>::max_children_per_face; ++i)
97 FETools::compute_face_embedding_matrices<dim, double>(*
this,
104 unsigned int target_row = 0;
105 for (
unsigned int d = 0; d < GeometryInfo<dim>::max_children_per_face; ++
d)
106 for (
unsigned int i = 0; i < face_embeddings[
d].
m(); ++i)
108 for (
unsigned int j = 0; j < face_embeddings[
d].
n(); ++j)
135 std::ostringstream namebuf;
136 namebuf <<
"FE_RaviartThomas<" << dim <<
">(" << this->degree - 1 <<
")";
138 return namebuf.str();
143std::unique_ptr<FiniteElement<dim, dim>>
146 return std::make_unique<FE_RaviartThomas<dim>>(*this);
160 const unsigned int n_interior_points = (deg > 0) ? cell_quadrature.
size() : 0;
165 const unsigned int face_no = 0;
167 unsigned int n_face_points = (dim > 1) ? 1 : 0;
169 for (
unsigned int d = 1;
d < dim; ++
d)
170 n_face_points *= deg + 1;
173 this->generalized_support_points.resize(
175 this->generalized_face_support_points[face_no].resize(n_face_points);
178 unsigned int current = 0;
182 QGauss<dim - 1> face_points(deg + 1);
186 boundary_weights.reinit(n_face_points,
legendre.n());
188 for (
unsigned int k = 0; k < n_face_points; ++k)
190 this->generalized_face_support_points[face_no][k] =
191 face_points.point(k);
195 for (
unsigned int i = 0; i <
legendre.n(); ++i)
197 boundary_weights(k, i) =
198 face_points.weight(k) *
199 legendre.compute_value(i, face_points.point(k));
206 for (; current < GeometryInfo<dim>::faces_per_cell * n_face_points;
211 this->generalized_support_points[current] = faces.
point(
222 std::unique_ptr<AnisotropicPolynomials<dim>> polynomials[dim];
223 for (
unsigned int dd = 0; dd < dim; ++dd)
225 std::vector<std::vector<Polynomials::Polynomial<double>>> poly(dim);
226 for (
unsigned int d = 0;
d < dim; ++
d)
230 polynomials[dd] = std::make_unique<AnisotropicPolynomials<dim>>(poly);
233 interior_weights.reinit(
236 for (
unsigned int k = 0; k < cell_quadrature.
size(); ++k)
238 this->generalized_support_points[current++] = cell_quadrature.
point(k);
239 for (
unsigned int i = 0; i < polynomials[0]->n(); ++i)
240 for (
unsigned int d = 0;
d < dim; ++
d)
241 interior_weights(k, i,
d) =
242 cell_quadrature.
weight(k) *
243 polynomials[
d]->compute_value(i, cell_quadrature.
point(k));
246 Assert(current == this->generalized_support_points.size(),
267 const unsigned int face_no = 0;
269 Assert(this->adjust_quad_dof_index_for_face_orientation_table[0]
270 .n_elements() == 8 * this->n_dofs_per_quad(face_no),
274 const unsigned int n = this->tensor_degree();
328 for (
unsigned int local_face_dof = 0;
329 local_face_dof < this->n_dofs_per_quad(face_no);
333 unsigned int i = local_face_dof % n;
334 unsigned int j = local_face_dof / n;
339 for (
unsigned int case_no = 0; case_no < 8; ++case_no)
346 if (((!face_orientation) && (!face_rotation)) ||
347 ((face_orientation) && (face_rotation)))
351 this->adjust_quad_dof_index_for_face_orientation_table[face_no](
352 local_face_dof, case_no) = j + i * n - local_face_dof;
357 this->adjust_quad_dof_index_for_face_orientation_table[face_no](
358 local_face_dof, case_no) = 0;
362 const unsigned int new_local_face_dof =
364 this->adjust_quad_dof_index_for_face_orientation_table[face_no](
365 local_face_dof, case_no);
367 i = new_local_face_dof % n;
368 j = new_local_face_dof / n;
375 if (!face_flip && face_rotation)
377 this->adjust_quad_dof_sign_for_face_orientation_table[face_no](
378 local_face_dof, case_no) = ((j % 2) == 1);
381 else if (face_flip && !face_rotation)
385 this->adjust_quad_dof_sign_for_face_orientation_table[face_no](
386 local_face_dof, case_no) = ((j % 2) == 1) != ((i % 2) == 1);
389 else if (face_flip && face_rotation)
391 this->adjust_quad_dof_sign_for_face_orientation_table[face_no](
392 local_face_dof, case_no) = ((i % 2) == 1);
401 if (!face_orientation)
402 this->adjust_quad_dof_sign_for_face_orientation_table[face_no](
403 local_face_dof, case_no) =
404 !this->adjust_quad_dof_sign_for_face_orientation_table[face_no](
405 local_face_dof, case_no);
418 for (
unsigned int i = 0; i < GeometryInfo<1>::max_children_per_cell; ++i)
419 this->restriction[0][i].
reinit(0, 0);
440 QGauss<dim - 1> q_base(this->degree);
441 const unsigned int n_face_points = q_base.size();
457 for (
unsigned int k = 0; k < q_face.
size(); ++k)
458 for (
unsigned int i = 0; i < this->n_dofs_per_cell(); ++i)
459 cached_values_on_face(i, k) = this->shape_value_component(
462 for (
unsigned int sub = 0; sub < GeometryInfo<dim>::max_children_per_face;
486 for (
unsigned int k = 0; k < n_face_points; ++k)
487 for (
unsigned int i_child = 0; i_child < this->n_dofs_per_cell();
489 for (
unsigned int i_face = 0;
490 i_face < this->n_dofs_per_face(face);
498 this->restriction[iso][child](
499 face * this->n_dofs_per_face(face) + i_face, i_child) +=
500 Utilities::fixed_power<dim - 1>(.5) * q_sub.
weight(k) *
501 cached_values_on_face(i_child, k) *
502 this->shape_value_component(
503 face * this->n_dofs_per_face(face) + i_face,
510 if (this->degree == 1)
515 std::unique_ptr<AnisotropicPolynomials<dim>> polynomials[dim];
516 for (
unsigned int dd = 0; dd < dim; ++dd)
518 std::vector<std::vector<Polynomials::Polynomial<double>>> poly(dim);
519 for (
unsigned int d = 0;
d < dim; ++
d)
525 polynomials[dd] = std::make_unique<AnisotropicPolynomials<dim>>(poly);
531 const unsigned int face_no = 0;
534 const unsigned int start_cell_dofs =
543 for (
unsigned int k = 0; k < q_cell.
size(); ++k)
544 for (
unsigned int i = 0; i < this->n_dofs_per_cell(); ++i)
545 for (
unsigned int d = 0;
d < dim; ++
d)
546 cached_values_on_cell(i, k,
d) =
547 this->shape_value_component(i, q_cell.
point(k),
d);
549 for (
unsigned int child = 0; child < GeometryInfo<dim>::max_children_per_cell;
557 for (
unsigned int k = 0; k < q_sub.
size(); ++k)
558 for (
unsigned int i_child = 0; i_child < this->n_dofs_per_cell();
560 for (
unsigned int d = 0;
d < dim; ++
d)
561 for (
unsigned int i_weight = 0; i_weight < polynomials[
d]->n();
564 this->restriction[iso][child](start_cell_dofs + i_weight * dim +
567 q_sub.
weight(k) * cached_values_on_cell(i_child, k,
d) *
568 polynomials[
d]->compute_value(i_weight, q_sub.
point(k));
576std::vector<unsigned int>
581 unsigned int dofs_per_face = 1;
582 for (
unsigned int d = 1;
d < dim; ++
d)
583 dofs_per_face *= deg + 1;
586 const unsigned int interior_dofs = dim * deg * dofs_per_face;
588 std::vector<unsigned int> dpo(dim + 1);
589 dpo[dim - 1] = dofs_per_face;
590 dpo[dim] = interior_dofs;
598std::pair<Table<2, bool>, std::vector<unsigned int>>
602 for (
unsigned int d = 0;
d < dim; ++
d)
603 for (
unsigned int i = 0; i < this->n_dofs_per_cell(); ++i)
604 constant_modes(
d, i) =
true;
605 std::vector<unsigned int> components;
606 for (
unsigned int d = 0;
d < dim; ++
d)
607 components.push_back(
d);
608 return std::pair<Table<2, bool>, std::vector<unsigned int>>(constant_modes,
622 const unsigned int face_index)
const
630 switch (this->degree)
643 return (face_index !=
665 std::vector<double> & nodal_values)
const
667 Assert(support_point_values.size() == this->generalized_support_points.size(),
669 this->generalized_support_points.size()));
670 Assert(nodal_values.size() == this->n_dofs_per_cell(),
672 Assert(support_point_values[0].size() == this->n_components(),
674 this->n_components()));
676 std::fill(nodal_values.begin(), nodal_values.end(), 0.);
678 const unsigned int n_face_points = boundary_weights.size(0);
680 for (
unsigned int k = 0; k < n_face_points; ++k)
681 for (
unsigned int i = 0; i < boundary_weights.size(1); ++i)
683 nodal_values[i + face * this->n_dofs_per_face(face)] +=
684 boundary_weights(k, i) *
685 support_point_values[face * n_face_points + k](
692 const unsigned int face_no = 0;
694 const unsigned int start_cell_dofs =
696 const unsigned int start_cell_points =
699 for (
unsigned int k = 0; k < interior_weights.size(0); ++k)
700 for (
unsigned int i = 0; i < interior_weights.size(1); ++i)
701 for (
unsigned int d = 0;
d < dim; ++
d)
702 nodal_values[start_cell_dofs + i * dim +
d] +=
703 interior_weights(k, i,
d) *
704 support_point_values[k + start_cell_points](
d);
720#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 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 Quadrature< dim > project_to_all_faces(const Quadrature< dim - 1 > &quadrature)
static Quadrature< dim > project_to_child(const Quadrature< dim > &quadrature, const unsigned int child_no)
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
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)
std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
void reference_cell(Triangulation< dim, spacedim > &tria, const ReferenceCell &reference_cell)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
bool get_bit(const unsigned char number, const unsigned int n)
void reinit(MatrixBlock< MatrixType > &v, const BlockSparsityPattern &p)
double legendre(unsigned int l, double x)
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)