From: Wolfgang Bangerth Date: Thu, 4 Feb 2021 00:24:23 +0000 (-0700) Subject: Get rid of the PolynomialType template argument of FE_Q_Base. X-Git-Tag: v9.3.0-rc1~513^2~1 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=07d7e9f0cad8dfee9df7a0a5f9ad3eeb8f70cf4d;p=dealii.git Get rid of the PolynomialType template argument of FE_Q_Base. --- diff --git a/include/deal.II/fe/fe_bernstein.h b/include/deal.II/fe/fe_bernstein.h index b734bc5d66..2064818ca9 100644 --- a/include/deal.II/fe/fe_bernstein.h +++ b/include/deal.II/fe/fe_bernstein.h @@ -60,8 +60,7 @@ DEAL_II_NAMESPACE_OPEN */ template -class FE_Bernstein - : public FE_Q_Base, dim, spacedim> +class FE_Bernstein : public FE_Q_Base { public: /** diff --git a/include/deal.II/fe/fe_q.h b/include/deal.II/fe/fe_q.h index 967c26a902..099faf2511 100644 --- a/include/deal.II/fe/fe_q.h +++ b/include/deal.II/fe/fe_q.h @@ -545,7 +545,7 @@ DEAL_II_NAMESPACE_OPEN * */ template -class FE_Q : public FE_Q_Base, dim, spacedim> +class FE_Q : public FE_Q_Base { public: /** diff --git a/include/deal.II/fe/fe_q_base.h b/include/deal.II/fe/fe_q_base.h index 1340d25ad4..b12f3fb380 100644 --- a/include/deal.II/fe/fe_q_base.h +++ b/include/deal.II/fe/fe_q_base.h @@ -34,18 +34,16 @@ DEAL_II_NAMESPACE_OPEN * functional as a stand-alone. The completion of definitions is left to the * derived classes. */ -template +template class FE_Q_Base : public FE_Poly { public: /** * Constructor. */ - FE_Q_Base(const PolynomialType & poly_space, - const FiniteElementData &fe_data, - const std::vector & restriction_is_additive_flags); + FE_Q_Base(const ScalarPolynomialsBase &poly_space, + const FiniteElementData & fe_data, + const std::vector & restriction_is_additive_flags); /** * Return the matrix interpolating from the given finite element to the @@ -330,7 +328,7 @@ protected: struct Implementation; // Declare implementation friend. - friend struct FE_Q_Base::Implementation; + friend struct FE_Q_Base::Implementation; private: /** diff --git a/include/deal.II/fe/fe_q_bubbles.h b/include/deal.II/fe/fe_q_bubbles.h index fd54bfd838..753d1086e1 100644 --- a/include/deal.II/fe/fe_q_bubbles.h +++ b/include/deal.II/fe/fe_q_bubbles.h @@ -88,8 +88,7 @@ DEAL_II_NAMESPACE_OPEN * the bubble enrichments in the middle of the cell. */ template -class FE_Q_Bubbles - : public FE_Q_Base, dim, spacedim> +class FE_Q_Bubbles : public FE_Q_Base { public: /** diff --git a/include/deal.II/fe/fe_q_dg0.h b/include/deal.II/fe/fe_q_dg0.h index d1ea8957b0..889d1e8ff8 100644 --- a/include/deal.II/fe/fe_q_dg0.h +++ b/include/deal.II/fe/fe_q_dg0.h @@ -235,8 +235,7 @@ DEAL_II_NAMESPACE_OPEN * */ template -class FE_Q_DG0 - : public FE_Q_Base, dim, spacedim> +class FE_Q_DG0 : public FE_Q_Base { public: /** diff --git a/include/deal.II/fe/fe_q_iso_q1.h b/include/deal.II/fe/fe_q_iso_q1.h index da212c964a..db8f9ba8d0 100644 --- a/include/deal.II/fe/fe_q_iso_q1.h +++ b/include/deal.II/fe/fe_q_iso_q1.h @@ -108,11 +108,7 @@ DEAL_II_NAMESPACE_OPEN * FE_Q_iso_Q1 with more than one subdivision does have less coupling. */ template -class FE_Q_iso_Q1 - : public FE_Q_Base< - TensorProductPolynomials>, - dim, - spacedim> +class FE_Q_iso_Q1 : public FE_Q_Base { public: /** diff --git a/source/fe/fe_bernstein.cc b/source/fe/fe_bernstein.cc index 635e92e65f..06f7b5bf6e 100644 --- a/source/fe/fe_bernstein.cc +++ b/source/fe/fe_bernstein.cc @@ -36,13 +36,13 @@ DEAL_II_NAMESPACE_OPEN template FE_Bernstein::FE_Bernstein(const unsigned int degree) - : FE_Q_Base, dim, spacedim>( - this->renumber_bases(degree), - FiniteElementData(this->get_dpo_vector(degree), - 1, - degree, - FiniteElementData::H1), - std::vector(1, false)) + : FE_Q_Base(this->renumber_bases(degree), + FiniteElementData(this->get_dpo_vector( + degree), + 1, + degree, + FiniteElementData::H1), + std::vector(1, false)) {} diff --git a/source/fe/fe_q.cc b/source/fe/fe_q.cc index 607e267e3b..9be753fa61 100644 --- a/source/fe/fe_q.cc +++ b/source/fe/fe_q.cc @@ -44,7 +44,7 @@ namespace internal return QGaussLobatto<1>(degree + 1).get_points(); else { - using FEQ = dealii::FE_Q_Base, 1, 1>; + using FEQ = dealii::FE_Q_Base<1, 1>; AssertThrow(false, FEQ::ExcFEQCannotHaveDegree0()); } return std::vector>(); @@ -57,7 +57,7 @@ namespace internal template FE_Q::FE_Q(const unsigned int degree) - : FE_Q_Base, dim, spacedim>( + : FE_Q_Base( TensorProductPolynomials( Polynomials::generate_complete_Lagrange_basis( internal::FE_Q::get_QGaussLobatto_points(degree))), @@ -74,7 +74,7 @@ FE_Q::FE_Q(const unsigned int degree) template FE_Q::FE_Q(const Quadrature<1> &points) - : FE_Q_Base, dim, spacedim>( + : FE_Q_Base( TensorProductPolynomials( Polynomials::generate_complete_Lagrange_basis(points.get_points())), FiniteElementData(this->get_dpo_vector(points.size() - 1), diff --git a/source/fe/fe_q_base.cc b/source/fe/fe_q_base.cc index 77a219bca4..25a4d0233f 100644 --- a/source/fe/fe_q_base.cc +++ b/source/fe/fe_q_base.cc @@ -82,8 +82,8 @@ namespace internal * A class with the same purpose as the similarly named class of the * Triangulation class. See there for more information. */ -template -struct FE_Q_Base::Implementation +template +struct FE_Q_Base::Implementation { /** * Initialize the hanging node constraints matrices. Called from the @@ -92,7 +92,7 @@ struct FE_Q_Base::Implementation template static void initialize_constraints(const std::vector> &, - FE_Q_Base &) + FE_Q_Base<1, spacedim> &) { // no constraints in 1d } @@ -101,13 +101,13 @@ struct FE_Q_Base::Implementation template static void initialize_constraints(const std::vector> & /*points*/, - FE_Q_Base &fe) + FE_Q_Base<2, spacedim> &fe) { const unsigned int dim = 2; unsigned int q_deg = fe.degree; - if (std::is_same>::value) + if (dynamic_cast *>( + &fe.get_poly_space()) != nullptr) q_deg = fe.degree - 1; // restricted to each face, the traces of the shape functions is an @@ -212,13 +212,13 @@ struct FE_Q_Base::Implementation template static void initialize_constraints(const std::vector> & /*points*/, - FE_Q_Base &fe) + FE_Q_Base<3, spacedim> &fe) { const unsigned int dim = 3; unsigned int q_deg = fe.degree; - if (std::is_same>::value) + if (dynamic_cast *>( + &fe.get_poly_space()) != nullptr) q_deg = fe.degree - 1; // For a detailed documentation of the interpolation see the @@ -410,28 +410,27 @@ struct FE_Q_Base::Implementation -template -FE_Q_Base::FE_Q_Base( - const PolynomialType & poly_space, - const FiniteElementData &fe_data, - const std::vector & restriction_is_additive_flags) +template +FE_Q_Base::FE_Q_Base( + const ScalarPolynomialsBase &poly_space, + const FiniteElementData & fe_data, + const std::vector & restriction_is_additive_flags) : FE_Poly( poly_space, fe_data, restriction_is_additive_flags, std::vector(1, std::vector(1, true))) - , q_degree(std::is_same>::value ? + , q_degree(dynamic_cast *>( + &poly_space) != nullptr ? this->degree - 1 : this->degree) {} -template +template void -FE_Q_Base::initialize( - const std::vector> &points) +FE_Q_Base::initialize(const std::vector> &points) { Assert(points[0][0] == 0, ExcMessage("The first support point has to be zero.")); @@ -505,16 +504,15 @@ FE_Q_Base::initialize( -template +template void -FE_Q_Base::get_interpolation_matrix( +FE_Q_Base::get_interpolation_matrix( const FiniteElement &x_source_fe, FullMatrix & interpolation_matrix) const { // go through the list of elements we can interpolate from - if (const FE_Q_Base *source_fe = - dynamic_cast *>( - &x_source_fe)) + if (const FE_Q_Base *source_fe = + dynamic_cast *>(&x_source_fe)) { // ok, source is a Q element, so we will be able to do the work Assert(interpolation_matrix.m() == this->n_dofs_per_cell(), @@ -604,9 +602,9 @@ FE_Q_Base::get_interpolation_matrix( -template +template void -FE_Q_Base::get_face_interpolation_matrix( +FE_Q_Base::get_face_interpolation_matrix( const FiniteElement &source_fe, FullMatrix & interpolation_matrix, const unsigned int face_no) const @@ -620,9 +618,9 @@ FE_Q_Base::get_face_interpolation_matrix( -template +template void -FE_Q_Base::get_subface_interpolation_matrix( +FE_Q_Base::get_subface_interpolation_matrix( const FiniteElement &x_source_fe, const unsigned int subface, FullMatrix & interpolation_matrix, @@ -633,9 +631,8 @@ FE_Q_Base::get_subface_interpolation_matrix( x_source_fe.n_dofs_per_face(face_no))); // see if source is a Q element - if (const FE_Q_Base *source_fe = - dynamic_cast *>( - &x_source_fe)) + if (const FE_Q_Base *source_fe = + dynamic_cast *>(&x_source_fe)) { // have this test in here since a table of size 2x0 reports its size as // 0x0 @@ -722,22 +719,21 @@ FE_Q_Base::get_subface_interpolation_matrix( -template +template bool -FE_Q_Base::hp_constraints_are_implemented() const +FE_Q_Base::hp_constraints_are_implemented() const { return true; } -template +template std::vector> -FE_Q_Base::hp_vertex_dof_identities( +FE_Q_Base::hp_vertex_dof_identities( const FiniteElement &fe_other) const { - if (dynamic_cast *>( - &fe_other) != nullptr) + if (dynamic_cast *>(&fe_other) != nullptr) { // there should be exactly one single DoF of each FE at a vertex, and they // should have identical value @@ -776,16 +772,15 @@ FE_Q_Base::hp_vertex_dof_identities( -template +template std::vector> -FE_Q_Base::hp_line_dof_identities( +FE_Q_Base::hp_line_dof_identities( const FiniteElement &fe_other) const { // we can presently only compute these identities if both FEs are FE_Qs or // if the other one is an FE_Nothing - if (const FE_Q_Base *fe_q_other = - dynamic_cast *>( - &fe_other)) + if (const FE_Q_Base *fe_q_other = + dynamic_cast *>(&fe_other)) { // dofs are located along lines, so two dofs are identical if they are // located at identical positions. if we had only equidistant points, we @@ -872,17 +867,16 @@ FE_Q_Base::hp_line_dof_identities( -template +template std::vector> -FE_Q_Base::hp_quad_dof_identities( +FE_Q_Base::hp_quad_dof_identities( const FiniteElement &fe_other, const unsigned int) const { // we can presently only compute these identities if both FEs are FE_Qs or // if the other one is an FE_Nothing - if (const FE_Q_Base *fe_q_other = - dynamic_cast *>( - &fe_other)) + if (const FE_Q_Base *fe_q_other = + dynamic_cast *>(&fe_other)) { // this works exactly like the line case above, except that now we have // to have two indices i1, i2 and j1, j2 to characterize the dofs on the @@ -949,9 +943,9 @@ FE_Q_Base::hp_quad_dof_identities( -template +template void -FE_Q_Base::initialize_unit_support_points( +FE_Q_Base::initialize_unit_support_points( const std::vector> &points) { const std::vector &index_map_inverse = @@ -974,9 +968,9 @@ FE_Q_Base::initialize_unit_support_points( -template +template void -FE_Q_Base::initialize_unit_face_support_points( +FE_Q_Base::initialize_unit_face_support_points( const std::vector> &points) { // no faces in 1d, so nothing to do @@ -1012,10 +1006,9 @@ FE_Q_Base::initialize_unit_face_support_points( -template +template void -FE_Q_Base:: - initialize_quad_dof_index_permutation() +FE_Q_Base::initialize_quad_dof_index_permutation() { // for 1D and 2D, do nothing if (dim < 3) @@ -1098,14 +1091,13 @@ FE_Q_Base:: -template +template unsigned int -FE_Q_Base::face_to_cell_index( - const unsigned int face_index, - const unsigned int face, - const bool face_orientation, - const bool face_flip, - const bool face_rotation) const +FE_Q_Base::face_to_cell_index(const unsigned int face_index, + const unsigned int face, + const bool face_orientation, + const bool face_flip, + const bool face_rotation) const { AssertIndexRange(face_index, this->n_dofs_per_face(face)); AssertIndexRange(face, GeometryInfo::faces_per_cell); @@ -1214,12 +1206,11 @@ FE_Q_Base::face_to_cell_index( -template +template std::vector -FE_Q_Base::get_dpo_vector( - const unsigned int degree) +FE_Q_Base::get_dpo_vector(const unsigned int degree) { - using FEQ = FE_Q_Base; + using FEQ = FE_Q_Base; AssertThrow(degree > 0, typename FEQ::ExcFEQCannotHaveDegree0()); std::vector dpo(dim + 1, 1U); for (unsigned int i = 1; i < dpo.size(); ++i) @@ -1229,9 +1220,9 @@ FE_Q_Base::get_dpo_vector( -template +template void -FE_Q_Base::initialize_constraints( +FE_Q_Base::initialize_constraints( const std::vector> &points) { Implementation::initialize_constraints(points, *this); @@ -1239,9 +1230,9 @@ FE_Q_Base::initialize_constraints( -template +template const FullMatrix & -FE_Q_Base::get_prolongation_matrix( +FE_Q_Base::get_prolongation_matrix( const unsigned int child, const RefinementCase &refinement_case) const { @@ -1438,9 +1429,9 @@ FE_Q_Base::get_prolongation_matrix( -template +template const FullMatrix & -FE_Q_Base::get_restriction_matrix( +FE_Q_Base::get_restriction_matrix( const unsigned int child, const RefinementCase &refinement_case) const { @@ -1582,9 +1573,9 @@ FE_Q_Base::get_restriction_matrix( //--------------------------------------------------------------------------- -template +template bool -FE_Q_Base::has_support_on_face( +FE_Q_Base::has_support_on_face( const unsigned int shape_index, const unsigned int face_index) const { @@ -1684,9 +1675,9 @@ FE_Q_Base::has_support_on_face( -template +template std::pair, std::vector> -FE_Q_Base::get_constant_modes() const +FE_Q_Base::get_constant_modes() const { Table<2, bool> constant_modes(1, this->n_dofs_per_cell()); // We here just care for the constant mode due to the polynomial space diff --git a/source/fe/fe_q_base.inst.in b/source/fe/fe_q_base.inst.in index 339f449ce1..4fb9aabcbc 100644 --- a/source/fe/fe_q_base.inst.in +++ b/source/fe/fe_q_base.inst.in @@ -18,19 +18,6 @@ for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS) { #if deal_II_dimension <= deal_II_space_dimension - template class FE_Q_Base, - deal_II_dimension, - deal_II_space_dimension>; - template class FE_Q_Base, - deal_II_dimension, - deal_II_space_dimension>; - template class FE_Q_Base, - deal_II_dimension, - deal_II_space_dimension>; - template class FE_Q_Base< - TensorProductPolynomials>, - deal_II_dimension, - deal_II_space_dimension>; + template class FE_Q_Base; #endif } diff --git a/source/fe/fe_q_bubbles.cc b/source/fe/fe_q_bubbles.cc index 8625b8e0ad..d449759f2e 100644 --- a/source/fe/fe_q_bubbles.cc +++ b/source/fe/fe_q_bubbles.cc @@ -186,15 +186,14 @@ namespace internal template FE_Q_Bubbles::FE_Q_Bubbles(const unsigned int q_degree) - : FE_Q_Base, dim, spacedim>( - TensorProductPolynomialsBubbles( - Polynomials::generate_complete_Lagrange_basis( - QGaussLobatto<1>(q_degree + 1).get_points())), - FiniteElementData(get_dpo_vector(q_degree), - 1, - q_degree + 1, - FiniteElementData::H1), - get_riaf_vector(q_degree)) + : FE_Q_Base(TensorProductPolynomialsBubbles( + Polynomials::generate_complete_Lagrange_basis( + QGaussLobatto<1>(q_degree + 1).get_points())), + FiniteElementData(get_dpo_vector(q_degree), + 1, + q_degree + 1, + FiniteElementData::H1), + get_riaf_vector(q_degree)) , n_bubbles((q_degree <= 1) ? 1 : dim) { Assert(q_degree > 0, @@ -226,7 +225,7 @@ FE_Q_Bubbles::FE_Q_Bubbles(const unsigned int q_degree) template FE_Q_Bubbles::FE_Q_Bubbles(const Quadrature<1> &points) - : FE_Q_Base, dim, spacedim>( + : FE_Q_Base( TensorProductPolynomialsBubbles( Polynomials::generate_complete_Lagrange_basis(points.get_points())), FiniteElementData(get_dpo_vector(points.size() - 1), @@ -455,8 +454,8 @@ FE_Q_Bubbles::has_support_on_face( if (shape_index >= this->n_dofs_per_cell() - n_bubbles) return false; else - return FE_Q_Base, dim, spacedim>:: - has_support_on_face(shape_index, face_index); + return FE_Q_Base::has_support_on_face(shape_index, + face_index); } diff --git a/source/fe/fe_q_dg0.cc b/source/fe/fe_q_dg0.cc index 647f841f4f..30a91675e9 100644 --- a/source/fe/fe_q_dg0.cc +++ b/source/fe/fe_q_dg0.cc @@ -33,15 +33,14 @@ DEAL_II_NAMESPACE_OPEN template FE_Q_DG0::FE_Q_DG0(const unsigned int degree) - : FE_Q_Base, dim, spacedim>( - TensorProductPolynomialsConst( - Polynomials::generate_complete_Lagrange_basis( - QGaussLobatto<1>(degree + 1).get_points())), - FiniteElementData(get_dpo_vector(degree), - 1, - degree, - FiniteElementData::L2), - get_riaf_vector(degree)) + : FE_Q_Base(TensorProductPolynomialsConst( + Polynomials::generate_complete_Lagrange_basis( + QGaussLobatto<1>(degree + 1).get_points())), + FiniteElementData(get_dpo_vector(degree), + 1, + degree, + FiniteElementData::L2), + get_riaf_vector(degree)) { Assert(degree > 0, ExcMessage("This element can only be used for polynomial degrees " @@ -61,7 +60,7 @@ FE_Q_DG0::FE_Q_DG0(const unsigned int degree) template FE_Q_DG0::FE_Q_DG0(const Quadrature<1> &points) - : FE_Q_Base, dim, spacedim>( + : FE_Q_Base( TensorProductPolynomialsConst( Polynomials::generate_complete_Lagrange_basis(points.get_points())), FiniteElementData(get_dpo_vector(points.size() - 1), @@ -226,8 +225,8 @@ FE_Q_DG0::get_interpolation_matrix( ExcDimensionMismatch(interpolation_matrix.m(), x_source_fe.n_dofs_per_cell())); - this->FE_Q_Base, dim, spacedim>:: - get_interpolation_matrix(x_source_fe, interpolation_matrix); + this->FE_Q_Base::get_interpolation_matrix( + x_source_fe, interpolation_matrix); } @@ -267,8 +266,8 @@ FE_Q_DG0::has_support_on_face( if (shape_index == this->n_dofs_per_cell() - 1) return true; else - return FE_Q_Base, dim, spacedim>:: - has_support_on_face(shape_index, face_index); + return FE_Q_Base::has_support_on_face(shape_index, + face_index); } diff --git a/source/fe/fe_q_iso_q1.cc b/source/fe/fe_q_iso_q1.cc index 23b99d0fe8..ed2fb26a31 100644 --- a/source/fe/fe_q_iso_q1.cc +++ b/source/fe/fe_q_iso_q1.cc @@ -32,10 +32,7 @@ DEAL_II_NAMESPACE_OPEN template FE_Q_iso_Q1::FE_Q_iso_Q1(const unsigned int subdivisions) - : FE_Q_Base< - TensorProductPolynomials>, - dim, - spacedim>( + : FE_Q_Base( TensorProductPolynomials>( Polynomials::generate_complete_Lagrange_basis_on_subdivisions( subdivisions,