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_raviart_thomas.h> 26 #include <deal.II/fe/fe_values.h> 27 #include <deal.II/fe/fe_tools.h> 32 DEAL_II_NAMESPACE_OPEN
44 std::vector<bool>(dim,true)))
46 Assert (dim >= 2, ExcImpossibleInDim(dim));
78 ref_case<RefinementCase<dim>::isotropic_refinement+1; ++ref_case)
82 for (
unsigned int i=0; i<nc; ++i)
83 this->
prolongation[ref_case-1][i].reinit (n_dofs, n_dofs);
89 for (
unsigned int i=0; i<GeometryInfo<dim>::max_children_per_face; ++i)
91 FETools::compute_face_embedding_matrices<dim,double>(*
this, face_embeddings, 0, 0);
94 unsigned int target_row=0;
95 for (
unsigned int d=0; d<GeometryInfo<dim>::max_children_per_face; ++d)
96 for (
unsigned int i=0; i<face_embeddings[d].m(); ++i)
98 for (
unsigned int j=0; j<face_embeddings[d].n(); ++j)
121 std::ostringstream namebuf;
122 namebuf <<
"FE_RaviartThomasNodal<" << dim <<
">(" << this->
degree-1 <<
")";
124 return namebuf.str();
150 unsigned int current = 0;
160 QGauss<dim-1> face_points (deg+1);
166 for (
unsigned int k=0;
170 ::DataSetDescriptor::face(0,
174 this->dofs_per_face));
176 current = this->dofs_per_face*GeometryInfo<dim>::faces_per_cell;
186 for (
unsigned int d=0; d<dim; ++d)
191 ((d==1) ? low : high));
193 ((d==1) ? low : high),
194 ((d==2) ? low : high));
195 Assert(dim<=3, ExcNotImplemented());
197 for (
unsigned int k=0; k<quadrature->
size(); ++k)
207 std::vector<unsigned int>
213 for (
unsigned int d=1; d<dim; ++d)
214 dofs_per_face *= deg+1;
220 std::vector<unsigned int> dpo(dim+1);
222 dpo[dim] = interior_dofs;
233 Assert (
false, ExcImpossibleInDim(1));
234 return std::vector<bool>();
245 for (
unsigned int d=2; d<dim; ++d)
246 dofs_per_face *= deg+1;
252 std::vector<bool> ret_val(dofs_per_cell,
false);
264 const unsigned int shape_index,
265 const unsigned int face_index)
const 268 ExcIndexRange (shape_index, 0, this->dofs_per_cell));
275 const unsigned int support_face = shape_index / this->
degree;
293 std::vector<double> &,
294 const std::vector<double> &)
const 296 Assert(
false, ExcNotImplemented());
303 std::vector<double> &local_dofs,
305 unsigned int offset)
const 310 ExcDimensionMismatch(local_dofs.size(),this->
dofs_per_cell));
312 ExcDimensionMismatch(values[0].size(),offset+this->
n_components()));
318 unsigned int fbase = 0;
320 for (; f<GeometryInfo<dim>::faces_per_cell;
331 const unsigned int istep = (this->
dofs_per_cell - fbase) / dim;
337 for (
unsigned int i=0; i<istep; ++i)
339 local_dofs[fbase+i] = values[fbase+i](offset+f);
344 Assert (fbase == this->dofs_per_cell, ExcInternalError());
351 std::vector<double> &local_dofs,
352 const VectorSlice<
const std::vector<std::vector<double> > > &values)
const 355 ExcDimensionMismatch(values.size(), this->
n_components()));
359 ExcDimensionMismatch(local_dofs.size(),this->
dofs_per_cell));
364 unsigned int fbase = 0;
366 for (; f<GeometryInfo<dim>::faces_per_cell;
376 const unsigned int istep = (this->
dofs_per_cell - fbase) / dim;
382 for (
unsigned int i=0; i<istep; ++i)
384 local_dofs[fbase+i] = values[f][fbase+i];
389 Assert (fbase == this->dofs_per_cell, ExcInternalError());
404 std::vector<std::pair<unsigned int, unsigned int> >
414 return std::vector<std::pair<unsigned int, unsigned int> > ();
417 Assert(
false, ExcNotImplemented());
418 return std::vector<std::pair<unsigned int, unsigned int> > ();
425 std::vector<std::pair<unsigned int, unsigned int> >
438 return std::vector<std::pair<unsigned int, unsigned int> >();
459 const unsigned int p = this->
degree-1;
460 const unsigned int q = fe_q_other->degree-1;
462 std::vector<std::pair<unsigned int, unsigned int> > identities;
465 for (
unsigned int i=0; i<p+1; ++i)
466 identities.push_back (std::make_pair(i,i));
468 else if (p%2==0 && q%2==0)
469 identities.push_back(std::make_pair(p/2,q/2));
475 Assert (
false, ExcNotImplemented());
476 return std::vector<std::pair<unsigned int, unsigned int> > ();
482 std::vector<std::pair<unsigned int, unsigned int> >
495 return std::vector<std::pair<unsigned int, unsigned int> >();
500 const unsigned int q = fe_q_other->dofs_per_quad;
502 std::vector<std::pair<unsigned int, unsigned int> > identities;
505 for (
unsigned int i=0; i<p; ++i)
506 identities.push_back (std::make_pair(i,i));
508 else if (p%2!=0 && q%2!=0)
509 identities.push_back(std::make_pair(p/2, q/2));
515 Assert (
false, ExcNotImplemented());
516 return std::vector<std::pair<unsigned int, unsigned int> > ();
529 if (this->degree < fe_q_other->
degree)
530 return FiniteElementDomination::this_element_dominates;
531 else if (this->degree == fe_q_other->degree)
532 return FiniteElementDomination::either_element_can_dominate;
534 return FiniteElementDomination::other_element_dominates;
537 Assert (
false, ExcNotImplemented());
538 return FiniteElementDomination::neither_element_dominates;
549 Assert (
false, ExcImpossibleInDim(1));
560 Assert (
false, ExcImpossibleInDim(1));
578 ExcInterpolationNotImplemented());
581 ExcDimensionMismatch (interpolation_matrix.
n(),
584 ExcDimensionMismatch (interpolation_matrix.
m(),
605 ExcInterpolationNotImplemented ());
613 Quadrature<dim-1> quad_face_support (source_fe.get_generalized_face_support_points ());
620 double eps = 2e-13*this->
degree*(dim-1);
630 const Point<dim> &p = face_projection.point (i);
644 if ( std::fabs(matrix_entry - 1.0) < eps )
646 if ( std::fabs(matrix_entry) < eps )
649 interpolation_matrix(i,j) = matrix_entry;
663 sum += interpolation_matrix(j,i);
675 const unsigned int subface,
685 ExcInterpolationNotImplemented());
688 ExcDimensionMismatch (interpolation_matrix.
n(),
691 ExcDimensionMismatch (interpolation_matrix.
m(),
712 ExcInterpolationNotImplemented ());
720 Quadrature<dim-1> quad_face_support (source_fe.get_generalized_face_support_points ());
727 double eps = 2e-13*this->
degree*(dim-1);
738 const Point<dim> &p = subface_projection.point (i);
751 if ( std::fabs(matrix_entry - 1.0) < eps )
753 if ( std::fabs(matrix_entry) < eps )
756 interpolation_matrix(i,j) = matrix_entry;
770 sum += interpolation_matrix(j,i);
780 #include "fe_raviart_thomas_nodal.inst" 783 DEAL_II_NAMESPACE_CLOSE
std::vector< Point< dim > > generalized_support_points
FullMatrix< double > interface_constraints
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
void invert(const FullMatrix< number2 > &M)
static unsigned int compute_n_pols(unsigned int degree)
#define AssertThrow(cond, exc)
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
static unsigned int n_children(const RefinementCase< dim > &refinement_case)
std::vector< std::vector< FullMatrix< double > > > prolongation
virtual std::string get_name() const
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)
void initialize_support_points(const unsigned int rt_degree)
virtual std::string get_name() const =0
virtual unsigned int face_to_cell_index(const unsigned int face_dof_index, const unsigned int face, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false) const
virtual bool hp_constraints_are_implemented() const
unsigned int size() const
const unsigned int dofs_per_cell
std::vector< Point< dim-1 > > generalized_face_support_points
static Quadrature< dim > project_to_all_faces(const SubQuadrature &quadrature)
unsigned int n_components() const
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
virtual FiniteElement< dim > * clone() const
virtual void interpolate(std::vector< double > &local_dofs, const std::vector< double > &values) const
virtual double shape_value_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
FullMatrix< double > inverse_node_matrix
static std::vector< bool > get_ria_vector(const unsigned int degree)