79 ref_case < RefinementCase<dim>::isotropic_refinement + 1;
82 const unsigned int nc =
85 for (
unsigned int i = 0; i < nc; ++i)
92 const unsigned int face_no = 0;
99 for (
unsigned int i = 0; i < GeometryInfo<dim>::max_children_per_face; ++i)
102 FETools::compute_face_embedding_matrices<dim, double>(*
this,
109 unsigned int target_row = 0;
110 for (
unsigned int d = 0; d < GeometryInfo<dim>::max_children_per_face; ++
d)
111 for (
unsigned int i = 0; i < face_embeddings[
d].
m(); ++i)
113 for (
unsigned int j = 0; j < face_embeddings[
d].
n(); ++j)
140 std::ostringstream namebuf;
141 namebuf <<
"FE_RaviartThomasNodal<" << dim <<
">(" << this->degree - 1 <<
")";
143 return namebuf.str();
148std::unique_ptr<FiniteElement<dim, dim>>
151 return std::make_unique<FE_RaviartThomasNodal<dim>>(*this);
168 const unsigned int face_no = 0;
170 this->generalized_support_points.resize(this->n_dofs_per_cell());
171 this->generalized_face_support_points[face_no].resize(
172 this->n_dofs_per_face(face_no));
175 unsigned int current = 0;
185 QGauss<dim - 1> face_points(deg + 1);
186 Assert(face_points.size() == this->n_dofs_per_face(face_no),
188 for (
unsigned int k = 0; k < this->n_dofs_per_face(face_no); ++k)
189 this->generalized_face_support_points[face_no][k] =
190 face_points.point(k);
194 for (
unsigned int k = 0; k < this->n_dofs_per_face(face_no) *
197 this->generalized_support_points[k] = faces.
point(
203 this->n_dofs_per_face(
218 for (
unsigned int d = 0;
d < dim; ++
d)
220 std::unique_ptr<QAnisotropic<dim>> quadrature;
224 quadrature = std::make_unique<QAnisotropic<dim>>(high);
228 std::make_unique<QAnisotropic<dim>>(((
d == 0) ? low : high),
229 ((
d == 1) ? low : high));
233 std::make_unique<QAnisotropic<dim>>(((
d == 0) ? low : high),
234 ((
d == 1) ? low : high),
235 ((
d == 2) ? low : high));
241 for (
unsigned int k = 0; k < quadrature->size(); ++k)
242 this->generalized_support_points[current++] = quadrature->point(k);
252 dim>::initialize_quad_dof_index_permutation_and_sign_change()
265std::vector<unsigned int>
270 unsigned int dofs_per_face = 1;
271 for (
unsigned int d = 1;
d < dim; ++
d)
272 dofs_per_face *= deg + 1;
275 const unsigned int interior_dofs = dim * deg * dofs_per_face;
277 std::vector<unsigned int> dpo(dim + 1);
278 dpo[dim - 1] = dofs_per_face;
279 dpo[dim] = interior_dofs;
291 return std::vector<bool>();
300 const unsigned int dofs_per_cell =
302 unsigned int dofs_per_face = deg + 1;
303 for (
unsigned int d = 2;
d < dim; ++
d)
304 dofs_per_face *= deg + 1;
310 std::vector<bool> ret_val(dofs_per_cell,
false);
323 const unsigned int shape_index,
324 const unsigned int face_index)
const
332 const unsigned int support_face = shape_index / this->degree;
353 std::vector<double> & nodal_values)
const
355 Assert(support_point_values.size() == this->generalized_support_points.size(),
357 this->generalized_support_points.size()));
358 Assert(nodal_values.size() == this->n_dofs_per_cell(),
360 Assert(support_point_values[0].size() == this->n_components(),
362 this->n_components()));
368 unsigned int fbase = 0;
370 for (; f < GeometryInfo<dim>::faces_per_cell;
371 ++f, fbase += this->n_dofs_per_face(f))
373 for (
unsigned int i = 0; i < this->n_dofs_per_face(f); ++i)
375 nodal_values[fbase + i] = support_point_values[fbase + i](
382 const unsigned int istep = (this->n_dofs_per_cell() - fbase) / dim;
386 while (fbase < this->n_dofs_per_cell())
388 for (
unsigned int i = 0; i < istep; ++i)
390 nodal_values[fbase + i] = support_point_values[fbase + i](f);
415std::vector<std::pair<unsigned int, unsigned int>>
425 return std::vector<std::pair<unsigned int, unsigned int>>();
427 return std::vector<std::pair<unsigned int, unsigned int>>();
431 return std::vector<std::pair<unsigned int, unsigned int>>();
438std::vector<std::pair<unsigned int, unsigned int>>
452 return std::vector<std::pair<unsigned int, unsigned int>>();
473 const unsigned int p = this->degree - 1;
474 const unsigned int q = fe_q_other->degree - 1;
476 std::vector<std::pair<unsigned int, unsigned int>> identities;
479 for (
unsigned int i = 0; i < p + 1; ++i)
480 identities.emplace_back(i, i);
482 else if (p % 2 == 0 && q % 2 == 0)
483 identities.emplace_back(p / 2, q / 2);
491 return std::vector<std::pair<unsigned int, unsigned int>>();
496 return std::vector<std::pair<unsigned int, unsigned int>>();
502std::vector<std::pair<unsigned int, unsigned int>>
505 const unsigned int face_no)
const
517 return std::vector<std::pair<unsigned int, unsigned int>>();
521 const unsigned int p = this->n_dofs_per_quad(face_no);
524 const unsigned int q = fe_q_other->n_dofs_per_quad(0);
526 std::vector<std::pair<unsigned int, unsigned int>> identities;
529 for (
unsigned int i = 0; i < p; ++i)
530 identities.emplace_back(i, i);
532 else if (p % 2 != 0 && q % 2 != 0)
533 identities.emplace_back(p / 2, q / 2);
541 return std::vector<std::pair<unsigned int, unsigned int>>();
546 return std::vector<std::pair<unsigned int, unsigned int>>();
555 const unsigned int codim)
const
565 if (this->degree < fe_rt_nodal_other->degree)
567 else if (this->degree == fe_rt_nodal_other->degree)
575 if (fe_nothing->is_dominating())
595 const unsigned int)
const
607 const unsigned int)
const
619 const unsigned int face_no)
const
626 &x_source_fe) !=
nullptr),
629 Assert(interpolation_matrix.
n() == this->n_dofs_per_face(face_no),
631 this->n_dofs_per_face(face_no)));
669 double eps = 2
e-13 * this->degree * (dim - 1);
681 const Point<dim> &p = face_projection.point(i);
683 for (
unsigned int j = 0; j < this->n_dofs_per_face(face_no); ++j)
685 double matrix_entry =
686 this->shape_value_component(this->face_to_cell_index(j, 0), p, 0);
699 interpolation_matrix(i, j) = matrix_entry;
712 for (
unsigned int i = 0; i < this->n_dofs_per_face(face_no); ++i)
713 sum += interpolation_matrix(j, i);
725 const unsigned int subface,
727 const unsigned int face_no)
const
734 &x_source_fe) !=
nullptr),
737 Assert(interpolation_matrix.
n() == this->n_dofs_per_face(face_no),
739 this->n_dofs_per_face(face_no)));
777 double eps = 2
e-13 * this->degree * (dim - 1);
791 const Point<dim> &p = subface_projection.point(i);
793 for (
unsigned int j = 0; j < this->n_dofs_per_face(face_no); ++j)
795 double matrix_entry =
796 this->shape_value_component(this->face_to_cell_index(j, 0), p, 0);
809 interpolation_matrix(i, j) = matrix_entry;
822 for (
unsigned int i = 0; i < this->n_dofs_per_face(face_no); ++i)
823 sum += interpolation_matrix(j, i);
833#include "fe_raviart_thomas_nodal.inst"
FullMatrix< double > inverse_node_matrix
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
void initialize_quad_dof_index_permutation_and_sign_change()
static std::vector< bool > get_ria_vector(const unsigned int degree)
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
void initialize_support_points(const unsigned int rt_degree)
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
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)
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
virtual std::string get_name() const =0
std::vector< std::vector< Point< dim - 1 > > > generalized_face_support_points
FullMatrix< double > interface_constraints
std::vector< std::vector< FullMatrix< double > > > prolongation
void invert(const FullMatrix< number2 > &M)
static unsigned int n_polynomials(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 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
SymmetricTensor< rank, dim, Number > sum(const SymmetricTensor< rank, dim, Number > &local, const MPI_Comm &mpi_communicator)
#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)
#define AssertThrow(cond, exc)
Expression fabs(const Expression &x)
std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
@ either_element_can_dominate
@ other_element_dominates
@ neither_element_dominates
void reference_cell(Triangulation< dim, spacedim > &tria, const ReferenceCell &reference_cell)
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
void reinit(MatrixBlock< MatrixType > &v, const BlockSparsityPattern &p)
static unsigned int n_children(const RefinementCase< dim > &refinement_case)