47 std::vector<unsigned int>
48 get_rt_dpo_vector(
const unsigned int dim,
const unsigned int degree)
50 std::vector<unsigned int> dpo(dim + 1);
53 unsigned int dofs_per_face = 1;
54 for (
unsigned int d = 1;
d < dim; ++
d)
55 dofs_per_face *= (degree + 1);
57 dpo[dim - 1] = dofs_per_face;
58 dpo[dim] = dim * degree * dofs_per_face;
91 this->n_dofs_per_cell());
93 const unsigned int face_no = 0;
100 for (
unsigned int i = 0; i < GeometryInfo<dim>::max_children_per_face; ++i)
103 FETools::compute_face_embedding_matrices<dim, double>(*
this,
110 unsigned int target_row = 0;
111 for (
unsigned int d = 0; d < GeometryInfo<dim>::max_children_per_face; ++d)
112 for (
unsigned int i = 0; i < face_embeddings[d].
m(); ++i)
114 for (
unsigned int j = 0; j < face_embeddings[d].
n(); ++j)
136 return "FE_RaviartThomasNodal<" + std::to_string(dim) +
">(" +
137 std::to_string(this->degree - 1) +
")";
142std::unique_ptr<FiniteElement<dim, dim>>
145 return std::make_unique<FE_RaviartThomasNodal<dim>>(*this);
158 dim>::initialize_quad_dof_index_permutation_and_sign_change()
164 const unsigned int n = this->degree;
165 const unsigned int face_no = 0;
167 for (
unsigned int local = 0; local < this->n_dofs_per_quad(face_no); ++local)
171 unsigned int i = local % n, j = local / n;
174 this->adjust_quad_dof_index_for_face_orientation_table[face_no](local,
178 this->adjust_quad_dof_index_for_face_orientation_table[face_no](local,
180 i + (n - 1 - j) * n - local;
182 this->adjust_quad_dof_index_for_face_orientation_table[face_no](local,
184 (n - 1 - j) + (n - 1 - i) * n - local;
186 this->adjust_quad_dof_index_for_face_orientation_table[face_no](local,
188 (n - 1 - i) + j * n - local;
190 this->adjust_quad_dof_index_for_face_orientation_table[face_no](local,
193 this->adjust_quad_dof_index_for_face_orientation_table[face_no](local,
195 j + (n - 1 - i) * n - local;
197 this->adjust_quad_dof_index_for_face_orientation_table[face_no](local,
199 (n - 1 - i) + (n - 1 - j) * n - local;
201 this->adjust_quad_dof_index_for_face_orientation_table[face_no](local,
203 (n - 1 - j) + i * n - local;
206 for (
unsigned int i = 0; i < 4; ++i)
207 this->adjust_quad_dof_sign_for_face_orientation_table[face_no](local,
217 const unsigned int shape_index,
218 const unsigned int face_index)
const
225 const unsigned int support_face = shape_index / this->n_dofs_per_face();
243 std::vector<double> & nodal_values)
const
245 Assert(support_point_values.size() == this->generalized_support_points.size(),
247 this->generalized_support_points.size()));
248 Assert(nodal_values.size() == this->n_dofs_per_cell(),
250 Assert(support_point_values[0].size() == this->n_components(),
252 this->n_components()));
256 unsigned int fbase = 0;
258 for (; f < GeometryInfo<dim>::faces_per_cell;
259 ++f, fbase += this->n_dofs_per_face(f))
261 for (
unsigned int i = 0; i < this->n_dofs_per_face(f); ++i)
263 nodal_values[fbase + i] = support_point_values[fbase + i](
269 const unsigned int istep = (this->n_dofs_per_cell() - fbase) / dim;
273 while (fbase < this->n_dofs_per_cell())
275 for (
unsigned int i = 0; i < istep; ++i)
277 nodal_values[fbase + i] = support_point_values[fbase + i](f);
302std::vector<std::pair<unsigned int, unsigned int>>
310 return std::vector<std::pair<unsigned int, unsigned int>>();
312 return std::vector<std::pair<unsigned int, unsigned int>>();
316 return std::vector<std::pair<unsigned int, unsigned int>>();
323std::vector<std::pair<unsigned int, unsigned int>>
334 return std::vector<std::pair<unsigned int, unsigned int>>();
349 const unsigned int p = this->degree - 1;
350 const unsigned int q = fe_q_other->degree - 1;
352 std::vector<std::pair<unsigned int, unsigned int>> identities;
355 for (
unsigned int i = 0; i < p + 1; ++i)
356 identities.emplace_back(i, i);
358 else if (p % 2 == 0 && q % 2 == 0)
359 identities.emplace_back(p / 2, q / 2);
367 return std::vector<std::pair<unsigned int, unsigned int>>();
372 return std::vector<std::pair<unsigned int, unsigned int>>();
378std::vector<std::pair<unsigned int, unsigned int>>
381 const unsigned int face_no)
const
390 return std::vector<std::pair<unsigned int, unsigned int>>();
393 const unsigned int p = this->n_dofs_per_quad(face_no);
396 const unsigned int q = fe_q_other->n_dofs_per_quad(0);
398 std::vector<std::pair<unsigned int, unsigned int>> identities;
401 for (
unsigned int i = 0; i < p; ++i)
402 identities.emplace_back(i, i);
404 else if (p % 2 != 0 && q % 2 != 0)
405 identities.emplace_back(p / 2, q / 2);
413 return std::vector<std::pair<unsigned int, unsigned int>>();
418 return std::vector<std::pair<unsigned int, unsigned int>>();
427 const unsigned int codim)
const
437 if (this->degree < fe_rt_nodal_other->degree)
439 else if (this->degree == fe_rt_nodal_other->degree)
447 if (fe_nothing->is_dominating())
467 const unsigned int)
const
479 const unsigned int)
const
491 const unsigned int face_no)
const
498 &x_source_fe) !=
nullptr),
501 Assert(interpolation_matrix.
n() == this->n_dofs_per_face(face_no),
503 this->n_dofs_per_face(face_no)));
529 double eps = 2e-13 * this->degree * (dim - 1);
540 const Point<dim> &p = face_projection.point(i);
542 for (
unsigned int j = 0; j < this->n_dofs_per_face(face_no); ++j)
544 double matrix_entry =
545 this->shape_value_component(this->face_to_cell_index(j, 0), p, 0);
550 if (std::fabs(matrix_entry - 1.0) < eps)
552 if (std::fabs(matrix_entry) < eps)
555 interpolation_matrix(i, j) = matrix_entry;
566 for (
unsigned int i = 0; i < this->n_dofs_per_face(face_no); ++i)
567 sum += interpolation_matrix(j, i);
569 Assert(std::fabs(sum - 1) < 2e-13 * this->degree * (dim - 1),
580 const unsigned int subface,
582 const unsigned int face_no)
const
588 &x_source_fe) !=
nullptr),
591 Assert(interpolation_matrix.
n() == this->n_dofs_per_face(face_no),
593 this->n_dofs_per_face(face_no)));
619 double eps = 2e-13 * this->degree * (dim - 1);
631 const Point<dim> &p = subface_projection.point(i);
633 for (
unsigned int j = 0; j < this->n_dofs_per_face(face_no); ++j)
635 double matrix_entry =
636 this->shape_value_component(this->face_to_cell_index(j, 0), p, 0);
641 if (std::fabs(matrix_entry - 1.0) < eps)
643 if (std::fabs(matrix_entry) < eps)
646 interpolation_matrix(i, j) = matrix_entry;
657 for (
unsigned int i = 0; i < this->n_dofs_per_face(face_no); ++i)
658 sum += interpolation_matrix(j, i);
660 Assert(std::fabs(sum - 1) < 2e-13 * this->degree * (dim - 1),
671 const unsigned int child,
678 "Prolongation matrices are only available for refined cells!"));
682 if (this->prolongation[refinement_case - 1][child].n() == 0)
684 std::lock_guard<std::mutex> lock(this->mutex);
687 if (this->prolongation[refinement_case - 1][child].n() ==
688 this->n_dofs_per_cell())
689 return this->prolongation[refinement_case - 1][child];
697 std::vector<std::vector<FullMatrix<double>>> isotropic_matrices(
699 isotropic_matrices.back().resize(
702 this->n_dofs_per_cell()));
705 isotropic_matrices.back());
721 return this->prolongation[refinement_case - 1][child];
729 const unsigned int child,
736 "Restriction matrices are only available for refined cells!"));
740 if (this->restriction[refinement_case - 1][child].n() == 0)
742 std::lock_guard<std::mutex> lock(this->mutex);
745 if (this->restriction[refinement_case - 1][child].n() ==
746 this->n_dofs_per_cell())
747 return this->restriction[refinement_case - 1][child];
755 std::vector<std::vector<FullMatrix<double>>> isotropic_matrices(
757 isotropic_matrices.back().resize(
760 this->n_dofs_per_cell()));
762 this_nonconst.
restriction[refinement_case - 1].swap(
763 isotropic_matrices.back());
779 return this->restriction[refinement_case - 1][child];
785#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
std::vector< Point< dim > > get_polynomial_support_points() const
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)
#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)