17 #include <deal.II/base/quadrature.h> 18 #include <deal.II/base/quadrature_lib.h> 19 #include <deal.II/lac/vector.h> 20 #include <deal.II/fe/fe.h> 21 #include <deal.II/fe/fe_dgq.h> 22 #include <deal.II/fe/fe_tools.h> 27 #include <deal.II/base/std_cxx14/memory.h> 30 DEAL_II_NAMESPACE_OPEN
39 std::vector<Point<1> >
40 get_QGaussLobatto_points (
const unsigned int degree)
45 return std::vector<Point<1> >(1,
Point<1>(0.5));
53 template <
int dim,
int spacedim>
77 template <
int dim,
int spacedim>
83 std::vector<bool>(
FiniteElementData<dim>(get_dpo_vector(polynomials.size()-1),1, polynomials.size()-1).dofs_per_cell, true),
96 template <
int dim,
int spacedim>
104 std::ostringstream namebuf;
107 <<
">(" << this->degree <<
")";
108 return namebuf.str();
113 template <
int dim,
int spacedim>
117 std::vector<double> &nodal_values)
const 120 this->get_unit_support_points().size());
122 nodal_values.size());
124 nodal_values.size());
126 for (
unsigned int i=0; i<this->dofs_per_cell; ++i)
130 nodal_values[i] = support_point_values[i](0);
136 template <
int dim,
int spacedim>
137 std::unique_ptr<FiniteElement<dim,spacedim> >
140 return std_cxx14::make_unique<FE_DGQ<dim, spacedim>>(*this);
149 template <
int dim,
int spacedim>
150 std::vector<unsigned int>
153 std::vector<unsigned int> dpo(dim+1, 0U);
155 for (
unsigned int i=1; i<dim; ++i)
162 template <
int dim,
int spacedim>
165 const char direction)
const 167 const unsigned int n = this->degree+1;
169 for (
unsigned int i=1; i<dim; ++i)
178 for (
unsigned int i=n; i>0;)
188 for (
unsigned int iz=0; iz<((dim>2) ? n:1); ++iz)
189 for (
unsigned int j=0; j<n; ++j)
190 for (
unsigned int i=0; i<n; ++i)
192 unsigned int k = n*i-j+n-1 + n*n*iz;
199 for (
unsigned int iz=0; iz<((dim>2) ? n:1); ++iz)
200 for (
unsigned int iy=0; iy<n; ++iy)
201 for (
unsigned int ix=0; ix<n; ++ix)
203 unsigned int k = n*ix-iy+n-1 + n*n*iz;
211 for (
unsigned int iz=0; iz<n; ++iz)
212 for (
unsigned int iy=0; iy<n; ++iy)
213 for (
unsigned int ix=0; ix<n; ++ix)
215 unsigned int k = n*(n*iy-iz+n-1) + ix;
223 for (
unsigned int iz=0; iz<n; ++iz)
224 for (
unsigned int iy=0; iy<n; ++iy)
225 for (
unsigned int ix=0; ix<n; ++ix)
227 unsigned int k = n*(n*iy-iz+n-1) + ix;
239 template <
int dim,
int spacedim>
250 typename FE::ExcInterpolationNotImplemented() );
257 Assert (interpolation_matrix.
m() == this->dofs_per_cell,
259 this->dofs_per_cell));
270 this->dofs_per_cell);
275 for (
unsigned int j=0; j<this->dofs_per_cell; ++j)
280 const Point<dim> p = this->unit_support_points[j];
281 for (
unsigned int i=0; i<this->dofs_per_cell; ++i)
282 cell_interpolation(j,i)
283 = this->poly_space.compute_value (i, p);
286 source_interpolation(j,i)
295 cell_interpolation.
mmult (interpolation_matrix,
296 source_interpolation);
299 for (
unsigned int i=0; i<this->dofs_per_cell; ++i)
301 if (std::fabs(interpolation_matrix(i,j)) < 1e-15)
302 interpolation_matrix(i,j) = 0.;
309 for (
unsigned int i=0; i<this->dofs_per_cell; ++i)
313 sum += interpolation_matrix(i,j);
315 Assert (std::fabs(
sum-1) < 5e-14*std::max(this->degree,1U)*dim,
322 template <
int dim,
int spacedim>
334 (void)interpolation_matrix;
337 typename FE::ExcInterpolationNotImplemented());
339 Assert (interpolation_matrix.
m() == 0,
342 Assert (interpolation_matrix.
n() == 0,
349 template <
int dim,
int spacedim>
362 (void)interpolation_matrix;
365 typename FE::ExcInterpolationNotImplemented());
367 Assert (interpolation_matrix.
m() == 0,
370 Assert (interpolation_matrix.
n() == 0,
377 template <
int dim,
int spacedim>
386 ExcMessage(
"Prolongation matrices are only available for refined cells!"));
391 if (this->prolongation[refinement_case-1][child].n() == 0)
396 if (this->prolongation[refinement_case-1][child].n() ==
398 return this->prolongation[refinement_case-1][child];
405 std::vector<std::vector<FullMatrix<double> > >
407 isotropic_matrices.back().
414 isotropic_matrices,
true);
415 this_nonconst.
prolongation[refinement_case-1].swap(isotropic_matrices.back());
438 return this->prolongation[refinement_case-1][child];
443 template <
int dim,
int spacedim>
452 ExcMessage(
"Restriction matrices are only available for refined cells!"));
457 if (this->restriction[refinement_case-1][child].n() == 0)
462 if (this->restriction[refinement_case-1][child].n() ==
464 return this->restriction[refinement_case-1][child];
471 std::vector<std::vector<FullMatrix<double> > >
473 isotropic_matrices.back().
480 isotropic_matrices,
true);
481 this_nonconst.
restriction[refinement_case-1].swap(isotropic_matrices.back());
504 return this->restriction[refinement_case-1][child];
509 template <
int dim,
int spacedim>
518 template <
int dim,
int spacedim>
519 std::vector<std::pair<unsigned int, unsigned int> >
528 std::vector<std::pair<unsigned int, unsigned int> > ();
533 template <
int dim,
int spacedim>
534 std::vector<std::pair<unsigned int, unsigned int> >
543 std::vector<std::pair<unsigned int, unsigned int> > ();
548 template <
int dim,
int spacedim>
549 std::vector<std::pair<unsigned int, unsigned int> >
558 std::vector<std::pair<unsigned int, unsigned int> > ();
563 template <
int dim,
int spacedim>
575 template <
int dim,
int spacedim>
578 const unsigned int face_index)
const 580 Assert (shape_index < this->dofs_per_cell,
585 unsigned int n = this->degree+1;
590 if (this->unit_support_points.empty())
596 bool support_points_on_boundary =
true;
597 for (
unsigned int d=0; d<dim; ++d)
598 if (std::abs(this->unit_support_points[0][d]) > 1e-13)
599 support_points_on_boundary =
false;
600 for (
unsigned int d=0; d<dim; ++d)
601 if (std::abs(this->unit_support_points.back()[d]-1.) > 1e-13)
602 support_points_on_boundary =
false;
603 if (support_points_on_boundary ==
false)
606 unsigned int n2 = n*n;
618 return (((shape_index == 0) && (face_index == 0)) ||
619 ((shape_index == this->degree) && (face_index == 1)));
624 if (face_index==0 && (shape_index % n) == 0)
626 if (face_index==1 && (shape_index % n) == this->degree)
628 if (face_index==2 && shape_index < n)
630 if (face_index==3 && shape_index >= this->dofs_per_cell-n)
637 const unsigned int in2 = shape_index % n2;
640 if (face_index==0 && (shape_index % n) == 0)
643 if (face_index==1 && (shape_index % n) == n-1)
646 if (face_index==2 && in2 < n )
649 if (face_index==3 && in2 >= n2-n)
652 if (face_index==4 && shape_index < n2)
655 if (face_index==5 && shape_index >= this->dofs_per_cell - n2)
668 template <
int dim,
int spacedim>
669 std::pair<Table<2,bool>, std::vector<unsigned int> >
673 constant_modes.
fill(
true);
674 return std::pair<Table<2,bool>, std::vector<unsigned int> >
675 (constant_modes, std::vector<unsigned int>(1, 0));
681 template <
int dim,
int spacedim>
693 template <
int dim,
int spacedim>
695 :
FE_DGQ<dim,spacedim>(
Polynomials::generate_complete_Lagrange_basis(points.get_points()))
704 template <
int dim,
int spacedim>
710 std::ostringstream namebuf;
711 bool equidistant =
true;
712 std::vector<double> points(this->degree+1);
714 std::vector<unsigned int> lexicographic = this->poly_space.get_numbering_inverse();
715 for (
unsigned int j=0; j<=this->degree; j++)
716 points[j] = this->unit_support_points[lexicographic[j]][0];
719 for (
unsigned int j=0; j<=this->degree; j++)
720 if (std::abs(points[j] - (
double)j/this->degree) > 1e-15)
725 if (this->degree == 0 && std::abs(points[0]-0.5) < 1e-15)
728 if (equidistant ==
true)
730 if (this->degree > 2)
731 namebuf <<
"FE_DGQArbitraryNodes<" <<
Utilities::dim_string(dim,spacedim) <<
">(QIterated(QTrapez()," << this->degree <<
"))";
734 return namebuf.str();
739 bool gauss_lobatto =
true;
740 for (
unsigned int j=0; j<=this->degree; j++)
741 if (points[j] != points_gl.
point(j)(0))
743 gauss_lobatto =
false;
747 if (gauss_lobatto ==
true)
750 return namebuf.str();
754 const QGauss<1> points_g(this->degree+1);
756 for (
unsigned int j=0; j<=this->degree; j++)
757 if (points[j] != points_g.
point(j)(0))
765 namebuf <<
"FE_DGQArbitraryNodes<" <<
Utilities::dim_string(dim,spacedim) <<
">(QGauss(" << this->degree+1 <<
"))";
766 return namebuf.str();
771 bool gauss_log =
true;
772 for (
unsigned int j=0; j<=this->degree; j++)
773 if (points[j] != points_glog.
point(j)(0))
779 if (gauss_log ==
true)
781 namebuf <<
"FE_DGQArbitraryNodes<" <<
Utilities::dim_string(dim,spacedim) <<
">(QGaussLog(" << this->degree+1 <<
"))";
782 return namebuf.str();
786 namebuf <<
"FE_DGQArbitraryNodes<" <<
Utilities::dim_string(dim,spacedim) <<
">(QUnknownNodes(" << this->degree+1 <<
"))";
787 return namebuf.str();
792 template <
int dim,
int spacedim>
796 std::vector<double> &nodal_values)
const 799 this->get_unit_support_points().size());
801 nodal_values.size());
803 nodal_values.size());
805 for (
unsigned int i=0; i<this->dofs_per_cell; ++i)
809 nodal_values[i] = support_point_values[i](0);
816 template <
int dim,
int spacedim>
817 std::unique_ptr<FiniteElement<dim,spacedim> >
821 std::vector<Point<1> > qpoints(this->degree+1);
822 std::vector<unsigned int> lexicographic = this->poly_space.get_numbering_inverse();
823 for (
unsigned int i=0; i<=this->degree; ++i)
824 qpoints[i] =
Point<1>(this->unit_support_points[lexicographic[i]][0]);
827 return std_cxx14::make_unique<FE_DGQArbitraryNodes<dim,spacedim>>(pquadrature);
834 template <
int dim,
int spacedim>
841 template <
int dim,
int spacedim>
842 std::pair<Table<2,bool>, std::vector<unsigned int> >
848 constant_modes(0,0) =
true;
849 return std::pair<Table<2,bool>, std::vector<unsigned int> >
850 (constant_modes, std::vector<unsigned int>(1, 0));
855 template <
int dim,
int spacedim>
865 template <
int dim,
int spacedim>
866 std::unique_ptr<FiniteElement<dim,spacedim> >
869 return std_cxx14::make_unique<FE_DGQLegendre<dim,spacedim>>(this->degree);
876 template <
int dim,
int spacedim>
878 :
FE_DGQ<dim,spacedim>(
Polynomials::HermiteLikeInterpolation::generate_complete_basis(degree))
883 template <
int dim,
int spacedim>
893 template <
int dim,
int spacedim>
894 std::unique_ptr<FiniteElement<dim,spacedim> >
897 return std_cxx14::make_unique<FE_DGQHermite<dim,spacedim>>(this->degree);
903 #include "fe_dgq.inst" 906 DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcFEHasNoSupportPoints()
virtual void convert_generalized_support_point_values_to_dof_values(const std::vector< Vector< double > > &support_point_values, std::vector< double > &nodal_values) const
virtual void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix) const
std::vector< std::vector< FullMatrix< double > > > restriction
#define AssertDimension(dim1, dim2)
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const
virtual FiniteElementDomination::Domination compare_for_face_domination(const FiniteElement< dim, spacedim > &fe_other) const
FE_DGQLegendre(const unsigned int degree)
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const
virtual void convert_generalized_support_point_values_to_dof_values(const std::vector< Vector< double > > &support_point_values, std::vector< double > &nodal_values) const
const unsigned int degree
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
const Point< dim > & point(const unsigned int i) const
#define AssertThrow(cond, exc)
virtual std::size_t memory_consumption() const
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
virtual const FullMatrix< double > & get_prolongation_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
std::vector< Point< dim > > unit_support_points
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const
static ::ExceptionBase & ExcMessage(std::string arg1)
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const
std::vector< std::vector< FullMatrix< double > > > prolongation
virtual std::string get_name() const
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
virtual const FullMatrix< double > & get_restriction_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
virtual void get_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
void rotate_indices(std::vector< unsigned int > &indices, const char direction) const
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
std::string dim_string(const int dim, const int spacedim)
SymmetricTensor< rank, dim, Number > sum(const SymmetricTensor< rank, dim, Number > &local, const MPI_Comm &mpi_communicator)
unsigned int size() const
const unsigned int dofs_per_cell
FE_DGQArbitraryNodes(const Quadrature< 1 > &points)
void mmult(FullMatrix< number2 > &C, const FullMatrix< number2 > &B, const bool adding=false) const
FE_DGQHermite(const unsigned int degree)
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_quad_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
void reinit_restriction_and_prolongation_matrices(const bool isotropic_restriction_only=false, const bool isotropic_prolongation_only=false)
static ::ExceptionBase & ExcNotImplemented()
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
double compute_value(const unsigned int i, const Point< dim > &p) const
TensorProductPolynomials< dim > poly_space
virtual std::string get_name() const
void fill(InputIterator entries, const bool C_style_indexing=true)
virtual bool hp_constraints_are_implemented() const
virtual std::string get_name() const
static ::ExceptionBase & ExcInternalError()
virtual std::string get_name() const