17 #include <deal.II/fe/fe_face.h> 18 #include <deal.II/fe/fe_poly_face.templates.h> 19 #include <deal.II/fe/fe_nothing.h> 20 #include <deal.II/fe/fe_tools.h> 21 #include <deal.II/base/quadrature_lib.h> 22 #include <deal.II/lac/householder.h> 25 DEAL_II_NAMESPACE_OPEN
30 std::vector<Point<1> >
31 get_QGaussLobatto_points (
const unsigned int degree)
36 return std::vector<Point<1> >(1,
Point<1>(0.5));
40 template <
int dim,
int spacedim>
46 std::vector<bool>(1,true))
49 const unsigned int codim = dim-1;
52 if (this->degree == 0)
53 for (
unsigned int d=0; d<codim; ++d)
57 std::vector<Point<1> > points = get_QGaussLobatto_points(
degree);
60 for (
unsigned int iz=0; iz <= ((codim>2) ? this->degree : 0) ; ++iz)
61 for (
unsigned int iy=0; iy <= ((codim>1) ? this->degree : 0) ; ++iy)
62 for (
unsigned int ix=0; ix<=this->
degree; ++ix)
82 for (
unsigned int i=0; i<n_face_dofs; ++i)
83 for (
unsigned int d=0; d<dim; ++d)
85 for (
unsigned int e=0, c=0; e<dim; ++e)
89 unsigned int renumber = i;
90 if (dim == 3 && d == 1)
104 template <
int dim,
int spacedim>
113 template <
int dim,
int spacedim>
120 std::ostringstream namebuf;
121 namebuf <<
"FE_FaceQ<" 123 <<
">(" << this->degree <<
")";
125 return namebuf.str();
130 template <
int dim,
int spacedim>
137 interpolation_matrix);
142 template <
int dim,
int spacedim>
146 const unsigned int subface,
151 Assert (interpolation_matrix.
n() == this->dofs_per_face,
153 this->dofs_per_face));
168 Assert (this->dofs_per_face <= source_fe->dofs_per_face,
170 ExcInterpolationNotImplemented ()));
173 const Quadrature<dim-1> face_quadrature (source_fe->get_unit_face_support_points ());
178 const double eps = 2e-13*(this->degree+1)*(dim-1);
182 for (
unsigned int i=0; i<source_fe->dofs_per_face; ++i)
184 const Point<dim-1> p =
187 face_quadrature.point(i)
192 for (
unsigned int j=0; j<this->dofs_per_face; ++j)
194 double matrix_entry = this->poly_space.compute_value (j, p);
199 if (std::fabs (matrix_entry - 1.0) < eps)
201 if (std::fabs (matrix_entry) < eps)
204 interpolation_matrix(i,j) = matrix_entry;
210 for (
unsigned int j=0; j<source_fe->dofs_per_face; ++j)
214 for (
unsigned int i=0; i<this->dofs_per_face; ++i)
215 sum += interpolation_matrix(j,i);
226 ExcInterpolationNotImplemented()));
231 template <
int dim,
int spacedim>
234 const unsigned int shape_index,
235 const unsigned int face_index)
const 237 return (face_index == (shape_index/this->dofs_per_face));
242 template <
int dim,
int spacedim>
243 std::vector<unsigned int>
246 std::vector<unsigned int> dpo(dim+1, 0U);
248 for (
unsigned int i=1; i<dim-1; ++i)
255 template <
int dim,
int spacedim>
264 template <
int dim,
int spacedim>
272 if (this->degree < fe_q_other->degree)
274 else if (this->degree == fe_q_other->degree)
281 if (fe_nothing->is_dominating())
299 template <
int dim,
int spacedim>
300 std::pair<Table<2,bool>, std::vector<unsigned int> >
304 for (
unsigned int i=0; i<this->dofs_per_cell; ++i)
305 constant_modes(0,i) =
true;
306 return std::pair<Table<2,bool>, std::vector<unsigned int> >
307 (constant_modes, std::vector<unsigned int>(1, 0));
314 template <
int spacedim>
318 std::vector<bool>(1,true),
331 template <
int spacedim>
340 template <
int spacedim>
347 std::ostringstream namebuf;
348 namebuf <<
"FE_FaceQ<" 350 <<
">(" << this->
degree <<
")";
352 return namebuf.str();
357 template <
int spacedim>
364 interpolation_matrix);
369 template <
int spacedim>
383 interpolation_matrix(0,0) = 1.;
388 template <
int spacedim>
391 const unsigned int shape_index,
392 const unsigned int face_index)
const 395 return (face_index == shape_index);
400 template <
int spacedim>
401 std::vector<unsigned int>
404 std::vector<unsigned int> dpo(2, 0U);
411 template <
int spacedim>
420 template <
int spacedim>
430 template <
int spacedim>
431 std::pair<Table<2,bool>, std::vector<unsigned int> >
436 constant_modes(0,i) =
true;
437 return std::pair<Table<2,bool>, std::vector<unsigned int> >
438 (constant_modes, std::vector<unsigned int>(1,0));
443 template <
int spacedim>
459 template <
int spacedim>
467 const ::internal::FEValues::MappingRelatedData<1, spacedim> &,
476 template <
int spacedim>
480 const unsigned int face,
484 const ::internal::FEValues::MappingRelatedData<1, spacedim> &,
488 const unsigned int foffset = face;
498 template <
int spacedim>
507 const ::internal::FEValues::MappingRelatedData<1, spacedim> &,
518 template <
int dim,
int spacedim>
524 std::vector<bool>(1,true))
529 template <
int dim,
int spacedim>
538 template <
int dim,
int spacedim>
545 std::ostringstream namebuf;
546 namebuf <<
"FE_FaceP<" 548 <<
">(" << this->degree <<
")";
550 return namebuf.str();
555 template <
int dim,
int spacedim>
558 const unsigned int shape_index,
559 const unsigned int face_index)
const 561 return (face_index == (shape_index/this->dofs_per_face));
566 template <
int dim,
int spacedim>
567 std::vector<unsigned int>
570 std::vector<unsigned int> dpo(dim+1, 0U);
572 for (
unsigned int i=1; i<dim-1; ++i)
574 dpo[dim-1] *= deg+1+i;
583 template <
int dim,
int spacedim>
592 template <
int dim,
int spacedim>
600 if (this->degree < fe_q_other->degree)
602 else if (this->degree == fe_q_other->degree)
609 if (fe_nothing->is_dominating())
628 template <
int dim,
int spacedim>
635 interpolation_matrix);
640 template <
int dim,
int spacedim>
644 const unsigned int subface,
649 Assert (interpolation_matrix.
n() == this->dofs_per_face,
651 this->dofs_per_face));
665 Assert (this->dofs_per_face <= source_fe->dofs_per_face,
667 ExcInterpolationNotImplemented ()));
672 const QGauss<dim-1> face_quadrature (source_fe->degree+1);
677 const double eps = 2e-13*(this->degree+1)*(dim-1);
681 for (
unsigned int k = 0; k < face_quadrature.size(); ++k)
683 const Point<dim-1> p =
685 face_quadrature.point(k) :
689 for (
unsigned int j = 0; j < source_fe->dofs_per_face; ++j)
690 mass (k , j) = source_fe->poly_space.compute_value(j, p);
700 for (
unsigned int i=0; i<this->dofs_per_face; ++i)
702 for (
unsigned int k = 0; k < face_quadrature.size(); ++k)
705 face_quadrature.point(k) :
708 v_in(k) = this->poly_space.compute_value(i, p);
714 for (
unsigned int j = 0; j < source_fe->dofs_per_face; ++j)
716 double matrix_entry = v_out(j);
721 if (std::fabs (matrix_entry - 1.0) < eps)
723 if (std::fabs (matrix_entry) < eps)
726 interpolation_matrix(j,i) = matrix_entry;
736 ExcInterpolationNotImplemented()));
741 template <
int dim,
int spacedim>
742 std::pair<Table<2,bool>, std::vector<unsigned int> >
746 for (
unsigned int face=0; face<GeometryInfo<dim>::faces_per_cell; ++face)
747 constant_modes(0, face*this->dofs_per_face) =
true;
748 return std::pair<Table<2,bool>, std::vector<unsigned int> >
749 (constant_modes, std::vector<unsigned int>(1, 0));
754 template <
int spacedim>
762 template <
int spacedim>
769 std::ostringstream namebuf;
770 namebuf <<
"FE_FaceP<" 772 <<
">(" << this->
degree <<
")";
774 return namebuf.str();
780 #include "fe_face.inst" 783 DEAL_II_NAMESPACE_CLOSE
static Point< dim > child_to_cell_coordinates(const Point< dim > &p, const unsigned int child_index, const RefinementCase< dim > refine_case=RefinementCase< dim >::isotropic_refinement)
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
Transformed quadrature weights.
static const unsigned int invalid_unsigned_int
static ::ExceptionBase & ExcLeastSquaresError(double arg1)
#define AssertDimension(dim1, dim2)
virtual FiniteElement< dim, spacedim > * clone() const
#define AssertIndexRange(index, range)
const unsigned int degree
virtual FiniteElement< dim, spacedim > * clone() const
virtual FiniteElementDomination::Domination compare_for_face_domination(const FiniteElement< dim, spacedim > &fe_other) const
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
#define AssertThrow(cond, exc)
virtual void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix) const
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const
std::vector< Point< dim > > unit_support_points
static ::ExceptionBase & ExcMessage(std::string arg1)
#define Assert(cond, exc)
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const
virtual std::string get_name() const
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
Abstract base class for mapping classes.
std::vector< Point< dim-1 > > unit_face_support_points
virtual bool hp_constraints_are_implemented() const
FE_FaceQ(const unsigned int p)
static std::vector< unsigned int > get_dpo_vector(const unsigned int deg)
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
Second derivatives of shape functions.
std::string dim_string(const int dim, const int spacedim)
virtual void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix) const
SymmetricTensor< rank, dim, Number > sum(const SymmetricTensor< rank, dim, Number > &local, const MPI_Comm &mpi_communicator)
const unsigned int dofs_per_cell
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
Shape function gradients.
static std::vector< unsigned int > get_dpo_vector(const unsigned int deg)
const unsigned int dofs_per_face
double least_squares(Vector< number2 > &dst, const Vector< number2 > &src) const
static ::ExceptionBase & ExcNotImplemented()
virtual std::string get_name() const
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const
Covariant transformation.
virtual FiniteElementDomination::Domination compare_for_face_domination(const FiniteElement< dim, spacedim > &fe_other) const
static ::ExceptionBase & ExcInternalError()
virtual bool hp_constraints_are_implemented() const