17 #include <deal.II/base/utilities.h> 18 #include <deal.II/base/quadrature.h> 19 #include <deal.II/base/quadrature_lib.h> 20 #include <deal.II/base/qprojector.h> 21 #include <deal.II/base/table.h> 22 #include <deal.II/grid/tria.h> 23 #include <deal.II/grid/tria_iterator.h> 24 #include <deal.II/dofs/dof_accessor.h> 25 #include <deal.II/fe/fe.h> 26 #include <deal.II/fe/mapping.h> 27 #include <deal.II/fe/fe_raviart_thomas.h> 28 #include <deal.II/fe/fe_values.h> 29 #include <deal.II/fe/fe_tools.h> 39 DEAL_II_NAMESPACE_OPEN
51 std::vector<bool>(dim,true)))
53 Assert (dim >= 2, ExcImpossibleInDim(dim));
92 for (
unsigned int i=0; i<GeometryInfo<dim>::max_children_per_face; ++i)
94 FETools::compute_face_embedding_matrices<dim,double>(*
this, face_embeddings, 0, 0);
97 unsigned int target_row=0;
98 for (
unsigned int d=0; d<GeometryInfo<dim>::max_children_per_face; ++d)
99 for (
unsigned int i=0; i<face_embeddings[d].m(); ++i)
101 for (
unsigned int j=0; j<face_embeddings[d].n(); ++j)
124 std::ostringstream namebuf;
125 namebuf <<
"FE_RaviartThomas<" << dim <<
">(" << this->
degree-1 <<
")";
127 return namebuf.str();
149 const unsigned int n_interior_points
150 = (deg>0) ? cell_quadrature.
size() : 0;
152 unsigned int n_face_points = (dim>1) ? 1 : 0;
154 for (
unsigned int d=1; d<dim; ++d)
155 n_face_points *= deg+1;
159 + n_interior_points);
163 unsigned int current = 0;
167 QGauss<dim-1> face_points (deg+1);
176 for (
unsigned int k=0; k<n_face_points; ++k)
182 for (
unsigned int i=0; i<legendre.n(); ++i)
185 = face_points.weight(k)
186 * legendre.compute_value(i, face_points.point(k));
191 for (; current<GeometryInfo<dim>::faces_per_cell*n_face_points;
204 std::vector<AnisotropicPolynomials<dim>* > polynomials(dim);
205 for (
unsigned int dd=0; dd<dim; ++dd)
207 std::vector<std::vector<Polynomials::Polynomial<double> > > poly(dim);
208 for (
unsigned int d=0; d<dim; ++d)
217 for (
unsigned int k=0; k<cell_quadrature.
size(); ++k)
220 for (
unsigned int i=0; i<polynomials[0]->n(); ++i)
221 for (
unsigned int d=0; d<dim; ++d)
223 * polynomials[d]->compute_value(i,cell_quadrature.
point(k));
226 for (
unsigned int d=0; d<dim; ++d)
227 delete polynomials[d];
242 for (
unsigned int i=0; i<GeometryInfo<1>::max_children_per_cell; ++i)
265 const unsigned int n_face_points = q_base.size();
268 for (
unsigned int face=0; face<GeometryInfo<dim>::faces_per_cell; ++face)
280 for (
unsigned int k=0; k<q_face.
size(); ++k)
286 for (
unsigned int sub=0; sub<GeometryInfo<dim>::max_children_per_face; ++sub)
294 const unsigned int child
310 for (
unsigned int k=0; k<n_face_points; ++k)
311 for (
unsigned int i_child = 0; i_child < this->
dofs_per_cell; ++i_child)
312 for (
unsigned int i_face = 0; i_face < this->
dofs_per_face; ++i_face)
319 this->
restriction[iso][child](face*this->dofs_per_face+i_face,
321 += Utilities::fixed_power<dim-1>(.5) * q_sub.
weight(k)
330 if (this->
degree == 1)
return;
335 std::vector<AnisotropicPolynomials<dim>* > polynomials(dim);
336 for (
unsigned int dd=0; dd<dim; ++dd)
338 std::vector<std::vector<Polynomials::Polynomial<double> > > poly(dim);
339 for (
unsigned int d=0; d<dim; ++d)
347 const unsigned int start_cell_dofs
354 for (
unsigned int k=0; k<q_cell.size(); ++k)
356 for (
unsigned int d=0; d<dim; ++d)
359 for (
unsigned int child=0; child<GeometryInfo<dim>::max_children_per_cell; ++child)
363 for (
unsigned int k=0; k<q_sub.
size(); ++k)
364 for (
unsigned int i_child = 0; i_child < this->
dofs_per_cell; ++i_child)
365 for (
unsigned int d=0; d<dim; ++d)
366 for (
unsigned int i_weight=0; i_weight<polynomials[d]->n(); ++i_weight)
368 this->
restriction[iso][child](start_cell_dofs+i_weight*dim+d,
372 * polynomials[d]->compute_value(i_weight, q_sub.
point(k));
376 for (
unsigned int d=0; d<dim; ++d)
377 delete polynomials[d];
383 std::vector<unsigned int>
389 for (
unsigned int d=1; d<dim; ++d)
390 dofs_per_face *= deg+1;
396 std::vector<unsigned int> dpo(dim+1);
398 dpo[dim] = interior_dofs;
406 std::pair<Table<2,bool>, std::vector<unsigned int> >
410 for (
unsigned int d=0; d<dim; ++d)
412 constant_modes(d,i) =
true;
414 for (
unsigned int d=0; d<dim; ++d)
415 components.push_back(d);
416 return std::pair<Table<2,bool>, std::vector<unsigned int> >
430 const unsigned int shape_index,
431 const unsigned int face_index)
const 434 ExcIndexRange (shape_index, 0, this->dofs_per_cell));
475 std::vector<double> &,
476 const std::vector<double> &)
const 478 Assert(
false, ExcNotImplemented());
485 std::vector<double> &local_dofs,
487 unsigned int offset)
const 492 ExcDimensionMismatch(local_dofs.size(),this->
dofs_per_cell));
494 ExcDimensionMismatch(values[0].size(),offset+this->
n_components()));
496 std::fill(local_dofs.begin(), local_dofs.end(), 0.);
499 for (
unsigned int face=0; face<GeometryInfo<dim>::faces_per_cell; ++face)
500 for (
unsigned int k=0; k<n_face_points; ++k)
512 for (
unsigned int d=0; d<dim; ++d)
513 local_dofs[start_cell_dofs+i*dim+d] +=
interior_weights(k,i,d) * values[k+start_cell_points](d+offset);
520 std::vector<double> &local_dofs,
521 const VectorSlice<
const std::vector<std::vector<double> > > &values)
const 524 ExcDimensionMismatch(values.size(), this->
n_components()));
528 ExcDimensionMismatch(local_dofs.size(),this->
dofs_per_cell));
530 std::fill(local_dofs.begin(), local_dofs.end(), 0.);
533 for (
unsigned int face=0; face<GeometryInfo<dim>::faces_per_cell; ++face)
534 for (
unsigned int k=0; k<n_face_points; ++k)
546 for (
unsigned int d=0; d<dim; ++d)
547 local_dofs[start_cell_dofs+i*dim+d] +=
interior_weights(k,i,d) * values[d][k+start_cell_points];
556 Assert (
false, ExcNotImplemented ());
563 #include "fe_raviart_thomas.inst" 566 DEAL_II_NAMESPACE_CLOSE
virtual std::string get_name() const
static Quadrature< dim > project_to_child(const Quadrature< dim > &quadrature, const unsigned int child_no)
virtual FiniteElement< dim > * clone() const
std::vector< std::vector< FullMatrix< double > > > restriction
const unsigned int components
std::vector< Point< dim > > generalized_support_points
FullMatrix< double > interface_constraints
static std::vector< Polynomial< double > > generate_complete_basis(const unsigned int degree)
const unsigned int degree
const Point< dim > & point(const unsigned int i) const
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const
void invert(const FullMatrix< number2 > &M)
friend class FE_RaviartThomas
std::vector< std::vector< FullMatrix< double > > > prolongation
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)
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::size_t memory_consumption() const
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
void initialize_restriction()
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)
std::vector< Tensor< 1, dim > > cached_values
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
void reinit_restriction_and_prolongation_matrices(const bool isotropic_restriction_only=false, const bool isotropic_prolongation_only=false)
virtual double shape_value_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
unsigned int size(const unsigned int i) const
virtual void interpolate(std::vector< double > &local_dofs, const std::vector< double > &values) const
FullMatrix< double > inverse_node_matrix
Table< 3, double > interior_weights
double weight(const unsigned int i) const
Table< 2, double > boundary_weights