17 #include <deal.II/base/quadrature_lib.h> 18 #include <deal.II/base/qprojector.h> 19 #include <deal.II/base/table.h> 20 #include <deal.II/grid/tria.h> 21 #include <deal.II/grid/tria_iterator.h> 22 #include <deal.II/dofs/dof_accessor.h> 23 #include <deal.II/fe/fe.h> 24 #include <deal.II/fe/mapping.h> 25 #include <deal.II/fe/fe_nothing.h> 26 #include <deal.II/fe/fe_raviart_thomas.h> 27 #include <deal.II/fe/fe_values.h> 28 #include <deal.II/fe/fe_tools.h> 31 #include <deal.II/base/std_cxx14/memory.h> 34 DEAL_II_NAMESPACE_OPEN
46 std::vector<bool>(dim,true)))
73 ref_case<RefinementCase<dim>::isotropic_refinement+1; ++ref_case)
77 for (
unsigned int i=0; i<nc; ++i)
78 this->
prolongation[ref_case-1][i].reinit (n_dofs, n_dofs);
84 for (
unsigned int i=0; i<GeometryInfo<dim>::max_children_per_face; ++i)
86 FETools::compute_face_embedding_matrices<dim,double>(*
this, face_embeddings, 0, 0);
89 unsigned int target_row=0;
90 for (
unsigned int d=0; d<GeometryInfo<dim>::max_children_per_face; ++d)
91 for (
unsigned int i=0; i<face_embeddings[d].m(); ++i)
93 for (
unsigned int j=0; j<face_embeddings[d].n(); ++j)
116 std::ostringstream namebuf;
117 namebuf <<
"FE_RaviartThomasNodal<" << dim <<
">(" << this->degree-1 <<
")";
119 return namebuf.str();
124 std::unique_ptr<FiniteElement<dim,dim> >
127 return std_cxx14::make_unique<FE_RaviartThomasNodal<dim>>(*this);
141 this->generalized_support_points.resize (this->dofs_per_cell);
142 this->generalized_face_support_points.resize (this->dofs_per_face);
145 unsigned int current = 0;
155 QGauss<dim-1> face_points (deg+1);
156 Assert (face_points.size() == this->dofs_per_face,
158 for (
unsigned int k=0; k<this->dofs_per_face; ++k)
159 this->generalized_face_support_points[k] = face_points.point(k);
161 for (
unsigned int k=0;
165 ::DataSetDescriptor::face(0,
169 this->dofs_per_face));
181 for (
unsigned int d=0; d<dim; ++d)
183 std::unique_ptr<QAnisotropic<dim>> quadrature;
187 quadrature = std_cxx14::make_unique<QAnisotropic<dim>>(high);
190 quadrature = std_cxx14::make_unique<QAnisotropic<dim>>(((d==0) ? low : high),
191 ((d==1) ? low : high));
194 quadrature = std_cxx14::make_unique<QAnisotropic<dim>>(((d==0) ? low : high),
195 ((d==1) ? low : high),
196 ((d==2) ? low : high));
202 for (
unsigned int k=0; k<quadrature->size(); ++k)
203 this->generalized_support_points[current++] = quadrature->point(k);
211 std::vector<unsigned int>
216 unsigned int dofs_per_face = 1;
217 for (
unsigned int d=1; d<dim; ++d)
218 dofs_per_face *= deg+1;
222 interior_dofs = dim*deg*dofs_per_face;
224 std::vector<unsigned int> dpo(dim+1);
225 dpo[dim-1] = dofs_per_face;
226 dpo[dim] = interior_dofs;
238 return std::vector<bool>();
248 unsigned int dofs_per_face = deg+1;
249 for (
unsigned int d=2; d<dim; ++d)
250 dofs_per_face *= deg+1;
256 std::vector<bool> ret_val(dofs_per_cell,
false);
258 i < dofs_per_cell; ++i)
268 const unsigned int shape_index,
269 const unsigned int face_index)
const 271 Assert (shape_index < this->dofs_per_cell,
279 const unsigned int support_face = shape_index / this->degree;
299 std::vector<double> &nodal_values)
const 301 Assert (support_point_values.size() == this->generalized_support_points.size(),
303 Assert (nodal_values.size() == this->dofs_per_cell,
312 unsigned int fbase = 0;
314 for (; f<GeometryInfo<dim>::faces_per_cell;
315 ++f, fbase+=this->dofs_per_face)
317 for (
unsigned int i=0; i<this->dofs_per_face; ++i)
325 const unsigned int istep = (this->dofs_per_cell - fbase) / dim;
329 while (fbase < this->dofs_per_cell)
331 for (
unsigned int i=0; i<istep; ++i)
333 nodal_values[fbase+i] = support_point_values[fbase+i](f);
354 std::vector<std::pair<unsigned int, unsigned int> >
364 return std::vector<std::pair<unsigned int, unsigned int> > ();
366 return std::vector<std::pair<unsigned int, unsigned int> > ();
370 return std::vector<std::pair<unsigned int, unsigned int> > ();
377 std::vector<std::pair<unsigned int, unsigned int> >
391 return std::vector<std::pair<unsigned int, unsigned int> >();
412 const unsigned int p = this->degree-1;
413 const unsigned int q = fe_q_other->
degree-1;
415 std::vector<std::pair<unsigned int, unsigned int> > identities;
418 for (
unsigned int i=0; i<p+1; ++i)
419 identities.emplace_back (i, i);
421 else if (p%2==0 && q%2==0)
422 identities.emplace_back(p/2, q/2);
430 return std::vector<std::pair<unsigned int, unsigned int> > ();
435 return std::vector<std::pair<unsigned int, unsigned int> > ();
441 std::vector<std::pair<unsigned int, unsigned int> >
455 return std::vector<std::pair<unsigned int, unsigned int> >();
459 const unsigned int p = this->dofs_per_quad;
462 std::vector<std::pair<unsigned int, unsigned int> > identities;
465 for (
unsigned int i=0; i<p; ++i)
466 identities.emplace_back (i, i);
468 else if (p%2!=0 && q%2!=0)
469 identities.emplace_back (p/2, q/2);
477 return std::vector<std::pair<unsigned int, unsigned int> > ();
482 return std::vector<std::pair<unsigned int, unsigned int> > ();
495 if (this->degree < fe_q_other->degree)
497 else if (this->degree == fe_q_other->degree)
505 if (fe_q_other->is_dominating())
559 ExcInterpolationNotImplemented());
561 Assert (interpolation_matrix.
n() == this->dofs_per_face,
563 this->dofs_per_face));
586 ExcInterpolationNotImplemented ());
601 double eps = 2
e-13*this->degree*(dim-1);
611 const Point<dim> &p = face_projection.point (i);
613 for (
unsigned int j=0; j<this->dofs_per_face; ++j)
616 = this->shape_value_component (this->face_to_cell_index(j, 0),
625 if ( std::fabs(matrix_entry - 1.0) < eps )
627 if ( std::fabs(matrix_entry) < eps )
630 interpolation_matrix(i,j) = matrix_entry;
643 for (
unsigned int i=0; i<this->dofs_per_face; ++i)
644 sum += interpolation_matrix(j,i);
646 Assert (std::fabs(
sum-1) < 2e-13*this->degree*(dim-1),
656 const unsigned int subface,
666 ExcInterpolationNotImplemented());
668 Assert (interpolation_matrix.
n() == this->dofs_per_face,
670 this->dofs_per_face));
693 ExcInterpolationNotImplemented ());
708 double eps = 2
e-13*this->degree*(dim-1);
719 const Point<dim> &p = subface_projection.point (i);
721 for (
unsigned int j=0; j<this->dofs_per_face; ++j)
724 = this->shape_value_component (this->face_to_cell_index(j, 0), p, 0);
732 if ( std::fabs(matrix_entry - 1.0) < eps )
734 if ( std::fabs(matrix_entry) < eps )
737 interpolation_matrix(i,j) = matrix_entry;
750 for (
unsigned int i=0; i<this->dofs_per_face; ++i)
751 sum += interpolation_matrix(j,i);
753 Assert (std::fabs(
sum-1) < 2e-13*this->degree*(dim-1),
761 #include "fe_raviart_thomas_nodal.inst" 764 DEAL_II_NAMESPACE_CLOSE
FullMatrix< double > interface_constraints
const std::vector< Point< dim-1 > > & get_generalized_face_support_points() const
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
virtual void convert_generalized_support_point_values_to_dof_values(const std::vector< Vector< double > > &support_point_values, std::vector< double > &nodal_values) const
FE_RaviartThomasNodal(const unsigned int p)
const unsigned int dofs_per_quad
const unsigned int degree
const Point< dim > & point(const unsigned int i) const
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
static unsigned int compute_n_pols(unsigned int degree)
#define AssertThrow(cond, exc)
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
std::vector< std::vector< FullMatrix< double > > > prolongation
virtual std::string get_name() const
void invert(const FullMatrix< number2 > &M)
void reinit(const TableIndices< N > &new_size, const bool omit_default_initialization=false)
#define Assert(cond, exc)
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 ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
void initialize_support_points(const unsigned int rt_degree)
virtual std::string get_name() const =0
virtual bool hp_constraints_are_implemented() const
SymmetricTensor< rank, dim, Number > sum(const SymmetricTensor< rank, dim, Number > &local, const MPI_Comm &mpi_communicator)
const unsigned int dofs_per_cell
static unsigned int n_children(const RefinementCase< dim > &refinement_case)
static Quadrature< dim > project_to_all_faces(const SubQuadrature &quadrature)
static void project_to_face(const SubQuadrature &quadrature, const unsigned int face_no, std::vector< Point< dim > > &q_points)
const unsigned int dofs_per_face
static ::ExceptionBase & ExcNotImplemented()
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
virtual std::unique_ptr< FiniteElement< dim, dim > > clone() const
FullMatrix< double > inverse_node_matrix
static ::ExceptionBase & ExcInternalError()
static std::vector< bool > get_ria_vector(const unsigned int degree)