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)))
85 for (
unsigned int i=0; i<GeometryInfo<dim>::max_children_per_face; ++i)
87 FETools::compute_face_embedding_matrices<dim,double>(*
this, face_embeddings, 0, 0);
90 unsigned int target_row=0;
91 for (
unsigned int d=0; d<GeometryInfo<dim>::max_children_per_face; ++d)
92 for (
unsigned int i=0; i<face_embeddings[d].m(); ++i)
94 for (
unsigned int j=0; j<face_embeddings[d].n(); ++j)
117 std::ostringstream namebuf;
118 namebuf <<
"FE_RaviartThomas<" << dim <<
">(" << this->degree-1 <<
")";
120 return namebuf.str();
142 const unsigned int n_interior_points
143 = (deg>0) ? cell_quadrature.
size() : 0;
145 unsigned int n_face_points = (dim>1) ? 1 : 0;
147 for (
unsigned int d=1; d<dim; ++d)
148 n_face_points *= deg+1;
152 + n_interior_points);
153 this->generalized_face_support_points.resize (n_face_points);
156 unsigned int current = 0;
160 QGauss<dim-1> face_points (deg+1);
164 boundary_weights.reinit(n_face_points, legendre.n());
169 for (
unsigned int k=0; k<n_face_points; ++k)
171 this->generalized_face_support_points[k] = face_points.point(k);
175 for (
unsigned int i=0; i<legendre.n(); ++i)
177 boundary_weights(k, i)
178 = face_points.weight(k)
179 * legendre.compute_value(i, face_points.point(k));
184 for (; current<GeometryInfo<dim>::faces_per_cell*n_face_points;
197 std::vector<AnisotropicPolynomials<dim>* > polynomials(dim);
198 for (
unsigned int dd=0; dd<dim; ++dd)
200 std::vector<std::vector<Polynomials::Polynomial<double> > > poly(dim);
201 for (
unsigned int d=0; d<dim; ++d)
208 interior_weights.reinit(
TableIndices<3>(n_interior_points, polynomials[0]->n(), dim));
210 for (
unsigned int k=0; k<cell_quadrature.
size(); ++k)
212 this->generalized_support_points[current++] = cell_quadrature.
point(k);
213 for (
unsigned int i=0; i<polynomials[0]->n(); ++i)
214 for (
unsigned int d=0; d<dim; ++d)
215 interior_weights(k,i,d) = cell_quadrature.
weight(k)
216 * polynomials[d]->compute_value(i,cell_quadrature.
point(k));
219 for (
unsigned int d=0; d<dim; ++d)
220 delete polynomials[d];
222 Assert (current == this->generalized_support_points.size(),
235 for (
unsigned int i=0; i<GeometryInfo<1>::max_children_per_cell; ++i)
236 this->restriction[0][i].reinit(0,0);
257 QGauss<dim-1> q_base (this->degree);
258 const unsigned int n_face_points = q_base.size();
261 for (
unsigned int face=0; face<GeometryInfo<dim>::faces_per_cell; ++face)
273 for (
unsigned int k=0; k<q_face.
size(); ++k)
274 for (
unsigned int i = 0; i < this->dofs_per_cell; ++i)
275 cached_values_on_face(i,k)
276 = this->shape_value_component(i, q_face.
point(k),
279 for (
unsigned int sub=0; sub<GeometryInfo<dim>::max_children_per_face; ++sub)
287 const unsigned int child
303 for (
unsigned int k=0; k<n_face_points; ++k)
304 for (
unsigned int i_child = 0; i_child < this->dofs_per_cell; ++i_child)
305 for (
unsigned int i_face = 0; i_face < this->dofs_per_face; ++i_face)
312 this->restriction[iso][child](face*this->dofs_per_face+i_face,
314 += Utilities::fixed_power<dim-1>(.5) * q_sub.
weight(k)
315 * cached_values_on_face(i_child, k)
316 * this->shape_value_component(face*this->dofs_per_face+i_face,
323 if (this->degree == 1)
return;
328 std::vector<AnisotropicPolynomials<dim>* > polynomials(dim);
329 for (
unsigned int dd=0; dd<dim; ++dd)
331 std::vector<std::vector<Polynomials::Polynomial<double> > > poly(dim);
332 for (
unsigned int d=0; d<dim; ++d)
340 const unsigned int start_cell_dofs
347 for (
unsigned int k=0; k<q_cell.
size(); ++k)
348 for (
unsigned int i = 0; i < this->dofs_per_cell; ++i)
349 for (
unsigned int d=0; d<dim; ++d)
350 cached_values_on_cell(i,k,d) = this->shape_value_component(i, q_cell.
point(k), d);
352 for (
unsigned int child=0; child<GeometryInfo<dim>::max_children_per_cell; ++child)
356 for (
unsigned int k=0; k<q_sub.
size(); ++k)
357 for (
unsigned int i_child = 0; i_child < this->dofs_per_cell; ++i_child)
358 for (
unsigned int d=0; d<dim; ++d)
359 for (
unsigned int i_weight=0; i_weight<polynomials[d]->n(); ++i_weight)
361 this->restriction[iso][child](start_cell_dofs+i_weight*dim+d,
364 * cached_values_on_cell(i_child, k, d)
365 * polynomials[d]->compute_value(i_weight, q_sub.
point(k));
369 for (
unsigned int d=0; d<dim; ++d)
370 delete polynomials[d];
376 std::vector<unsigned int>
381 unsigned int dofs_per_face = 1;
382 for (
unsigned int d=1; d<dim; ++d)
383 dofs_per_face *= deg+1;
387 interior_dofs = dim*deg*dofs_per_face;
389 std::vector<unsigned int> dpo(dim+1);
390 dpo[dim-1] = dofs_per_face;
391 dpo[dim] = interior_dofs;
399 std::pair<Table<2,bool>, std::vector<unsigned int> >
403 for (
unsigned int d=0; d<dim; ++d)
404 for (
unsigned int i=0; i<this->dofs_per_cell; ++i)
405 constant_modes(d,i) =
true;
406 std::vector<unsigned int> components;
407 for (
unsigned int d=0; d<dim; ++d)
408 components.push_back(d);
409 return std::pair<Table<2,bool>, std::vector<unsigned int> >
410 (constant_modes, components);
423 const unsigned int shape_index,
424 const unsigned int face_index)
const 426 Assert (shape_index < this->dofs_per_cell,
434 switch (this->degree)
468 std::vector<double> &nodal_values)
const 470 Assert (support_point_values.size() == this->generalized_support_points.size(),
472 Assert (nodal_values.size() == this->dofs_per_cell,
477 std::fill(nodal_values.begin(), nodal_values.end(), 0.);
479 const unsigned int n_face_points = boundary_weights.size(0);
480 for (
unsigned int face=0; face<GeometryInfo<dim>::faces_per_cell; ++face)
481 for (
unsigned int k=0; k<n_face_points; ++k)
482 for (
unsigned int i=0; i<boundary_weights.size(1); ++i)
484 nodal_values[i+face*this->dofs_per_face] += boundary_weights(k,i)
491 for (
unsigned int k=0; k<interior_weights.size(0); ++k)
492 for (
unsigned int i=0; i<interior_weights.size(1); ++i)
493 for (
unsigned int d=0; d<dim; ++d)
494 nodal_values[start_cell_dofs+i*dim+d] += interior_weights(k,i,d) * support_point_values[k+start_cell_points](d);
504 std::vector<double> &,
505 const std::vector<double> &)
const 514 std::vector<double> &local_dofs,
516 const unsigned int offset)
const 518 Assert (values.size() == this->generalized_support_points.size(),
520 Assert (local_dofs.size() == this->dofs_per_cell,
525 std::fill(local_dofs.begin(), local_dofs.end(), 0.);
527 const unsigned int n_face_points = boundary_weights.size(0);
528 for (
unsigned int face=0; face<GeometryInfo<dim>::faces_per_cell; ++face)
529 for (
unsigned int k=0; k<n_face_points; ++k)
530 for (
unsigned int i=0; i<boundary_weights.size(1); ++i)
532 local_dofs[i+face*this->dofs_per_face] += boundary_weights(k,i)
539 for (
unsigned int k=0; k<interior_weights.size(0); ++k)
540 for (
unsigned int i=0; i<interior_weights.size(1); ++i)
541 for (
unsigned int d=0; d<dim; ++d)
542 local_dofs[start_cell_dofs+i*dim+d] += interior_weights(k,i,d) * values[k+start_cell_points](d+offset);
549 std::vector<double> &local_dofs,
550 const VectorSlice<
const std::vector<std::vector<double> > > &values)
const 554 Assert (values[0].size() == this->generalized_support_points.size(),
556 Assert (local_dofs.size() == this->dofs_per_cell,
559 std::fill(local_dofs.begin(), local_dofs.end(), 0.);
561 const unsigned int n_face_points = boundary_weights.size(0);
562 for (
unsigned int face=0; face<GeometryInfo<dim>::faces_per_cell; ++face)
563 for (
unsigned int k=0; k<n_face_points; ++k)
564 for (
unsigned int i=0; i<boundary_weights.size(1); ++i)
566 local_dofs[i+face*this->dofs_per_face] += boundary_weights(k,i)
573 for (
unsigned int k=0; k<interior_weights.size(0); ++k)
574 for (
unsigned int i=0; i<interior_weights.size(1); ++i)
575 for (
unsigned int d=0; d<dim; ++d)
576 local_dofs[start_cell_dofs+i*dim+d] += interior_weights(k,i,d) * values[d][k+start_cell_points];
592 #include "fe_raviart_thomas.inst" 595 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
virtual void interpolate(std::vector< double > &local_dofs, const std::vector< double > &values) const 1
FullMatrix< double > interface_constraints
static std::vector< Polynomial< double > > generate_complete_basis(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
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
friend class FE_RaviartThomas
virtual void convert_generalized_support_point_values_to_nodal_values(const std::vector< Vector< double > > &support_point_values, std::vector< double > &nodal_values) const
std::vector< std::vector< FullMatrix< double > > > prolongation
void invert(const FullMatrix< number2 > &M)
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)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
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
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
void reinit_restriction_and_prolongation_matrices(const bool isotropic_restriction_only=false, const bool isotropic_prolongation_only=false)
static ::ExceptionBase & ExcNotImplemented()
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
FullMatrix< double > inverse_node_matrix
double weight(const unsigned int i) const
static ::ExceptionBase & ExcInternalError()