]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Make get_gauss_type_quadrature() and get_nodal_type_quadrature() members of Reference... 11664/head
authorWolfgang Bangerth <bangerth@colostate.edu>
Wed, 3 Feb 2021 04:08:41 +0000 (21:08 -0700)
committerWolfgang Bangerth <bangerth@colostate.edu>
Wed, 3 Feb 2021 04:33:01 +0000 (21:33 -0700)
include/deal.II/fe/fe_tools.templates.h
include/deal.II/grid/reference_cell.h
source/fe/mapping_fe_field.cc
source/grid/reference_cell.cc
source/grid/reference_cell.inst.in
tests/simplex/fe_p_bubbles_02.cc

index b2142de6293a7e0594939afa176470fcb7523df3..efdea4bd2f2019e585f8c46b57f4e006e2c3d906 100644 (file)
@@ -1598,7 +1598,7 @@ namespace FETools
     const unsigned int n1 = fe1.n_dofs_per_cell();
     const unsigned int n2 = fe2.n_dofs_per_cell();
 
-    const auto reference_cell_type = fe1.reference_cell_type();
+    const ReferenceCell::Type reference_cell_type = fe1.reference_cell_type();
 
     Assert(fe1.reference_cell_type() == fe2.reference_cell_type(),
            ExcNotImplemented());
@@ -1615,8 +1615,7 @@ namespace FETools
       std::max(fe1.tensor_degree(), fe2.tensor_degree());
     Assert(degree != numbers::invalid_unsigned_int, ExcNotImplemented());
     const auto quadrature =
-      ReferenceCell::get_gauss_type_quadrature<dim>(reference_cell_type,
-                                                    degree + 1);
+      reference_cell_type.get_gauss_type_quadrature<dim>(degree + 1);
 
     // Set up FEValues.
     const UpdateFlags flags =
@@ -1809,7 +1808,8 @@ namespace FETools
                    ExcDimensionMismatch(matrices[i].m(), n));
           }
 
-        const auto reference_cell_type = fe.reference_cell_type();
+        const ReferenceCell::Type reference_cell_type =
+          fe.reference_cell_type();
 
         // Set up meshes, one with a single
         // reference cell and refine it once
@@ -1824,8 +1824,7 @@ namespace FETools
           reference_cell_type
             .template get_default_linear_mapping<dim, spacedim>();
         const auto &q_fine =
-          ReferenceCell::get_gauss_type_quadrature<dim>(reference_cell_type,
-                                                        degree + 1);
+          reference_cell_type.get_gauss_type_quadrature<dim>(degree + 1);
         const unsigned int nq = q_fine.size();
 
         FEValues<dim, spacedim> fine(mapping,
index 892fec8364f245441de5235c3a55f09673f9f8c5..84b587b9b32e681f7e4bd799f3518ca7b9078e4c 100644 (file)
@@ -208,6 +208,28 @@ namespace ReferenceCell
     const Mapping<dim, spacedim> &
     get_default_linear_mapping() const;
 
+    /**
+     * Return a Gauss-type quadrature matching the given reference cell (QGauss,
+     * Simplex::QGauss, Simplex::QGaussPyramid, Simplex::QGaussWedge).
+     *
+     * @param[in] n_points_1D The number of quadrature points in each direction
+     * (QGauss) or an indication of what polynomial degree needs to be
+     * integrated exactly for the other types.
+     */
+    template <int dim>
+    Quadrature<dim>
+    get_gauss_type_quadrature(const unsigned n_points_1D) const;
+
+    /**
+     * Return a quadrature rule with the support points of the given reference
+     * cell.
+     *
+     * @note The weights of the quadrature object are left unfilled.
+     */
+    template <int dim>
+    const Quadrature<dim> &
+    get_nodal_type_quadrature() const;
+
     /**
      * Return a text representation of the reference cell represented by the
      * current object.
@@ -728,27 +750,6 @@ namespace ReferenceCell
   const Mapping<dim, spacedim> &
   get_default_linear_mapping(const Triangulation<dim, spacedim> &triangulation);
 
-  /**
-   * Return a Gauss-type quadrature matching the given reference cell(QGauss,
-   * Simplex::QGauss, Simplex::QGaussPyramid, Simplex::QGaussWedge) and
-   * @p n_points_1D the number of quadrature points in each direction (QGuass)
-   * or the indication of what polynomial degree to be integrated exactly.
-   */
-  template <int dim>
-  Quadrature<dim>
-  get_gauss_type_quadrature(const Type &   reference_cell,
-                            const unsigned n_points_1D);
-
-  /**
-   * Return a quadrature rule with the support points of the given reference
-   * cell.
-   *
-   * @note The weights are not filled.
-   */
-  template <int dim>
-  const Quadrature<dim> &
-  get_nodal_type_quadrature(const Type &reference_cell);
-
   namespace internal
   {
     /**
index 3b9712365a63d5e8035756891061004b0b768cea..91e4c9ea801d2309ff71e777219698b9d5cb3d19 100644 (file)
@@ -207,8 +207,9 @@ MappingFEField<dim, spacedim, VectorType, void>::MappingFEField(
                 true))
   , fe_to_real(fe_mask.size(), numbers::invalid_unsigned_int)
   , fe_values(this->euler_dof_handler->get_fe(),
-              ReferenceCell::get_nodal_type_quadrature<dim>(
-                this->euler_dof_handler->get_fe().reference_cell_type()),
+              this->euler_dof_handler->get_fe()
+                .reference_cell_type()
+                .template get_nodal_type_quadrature<dim>(),
               update_values)
 {
   unsigned int size = 0;
@@ -236,8 +237,9 @@ MappingFEField<dim, spacedim, VectorType, void>::MappingFEField(
                 true))
   , fe_to_real(fe_mask.size(), numbers::invalid_unsigned_int)
   , fe_values(this->euler_dof_handler->get_fe(),
-              ReferenceCell::get_nodal_type_quadrature<dim>(
-                this->euler_dof_handler->get_fe().reference_cell_type()),
+              this->euler_dof_handler->get_fe()
+                .reference_cell_type()
+                .template get_nodal_type_quadrature<dim>(),
               update_values)
 {
   unsigned int size = 0;
@@ -276,8 +278,9 @@ MappingFEField<dim, spacedim, VectorType, void>::MappingFEField(
                 true))
   , fe_to_real(fe_mask.size(), numbers::invalid_unsigned_int)
   , fe_values(this->euler_dof_handler->get_fe(),
-              ReferenceCell::get_nodal_type_quadrature<dim>(
-                this->euler_dof_handler->get_fe().reference_cell_type()),
+              this->euler_dof_handler->get_fe()
+                .reference_cell_type()
+                .template get_nodal_type_quadrature<dim>(),
               update_values)
 {
   unsigned int size = 0;
@@ -313,8 +316,9 @@ MappingFEField<dim, spacedim, VectorType, void>::MappingFEField(
   , fe_mask(mapping.fe_mask)
   , fe_to_real(mapping.fe_to_real)
   , fe_values(mapping.euler_dof_handler->get_fe(),
-              ReferenceCell::get_nodal_type_quadrature<dim>(
-                this->euler_dof_handler->get_fe().reference_cell_type()),
+              this->euler_dof_handler->get_fe()
+                .reference_cell_type()
+                .template get_nodal_type_quadrature<dim>(),
               update_values)
 {}
 
index c187cdd4ea278bf752fab9f530d5e15e5b8b50a8..ef0d1afcdc4146ed4ff6476733059145a792ae87 100644 (file)
@@ -179,18 +179,17 @@ namespace ReferenceCell
 
   template <int dim>
   Quadrature<dim>
-  get_gauss_type_quadrature(const Type &   reference_cell,
-                            const unsigned n_points_1D)
+  Type::get_gauss_type_quadrature(const unsigned n_points_1D) const
   {
-    AssertDimension(dim, reference_cell.get_dimension());
+    AssertDimension(dim, get_dimension());
 
-    if (reference_cell == Type::get_hypercube<dim>())
+    if (is_hyper_cube())
       return QGauss<dim>(n_points_1D);
-    else if (reference_cell == Type::Tri || reference_cell == Type::Tet)
+    else if (is_simplex())
       return Simplex::QGauss<dim>(n_points_1D);
-    else if (reference_cell == Type::Pyramid)
+    else if (*this == Type::Pyramid)
       return Simplex::QGaussPyramid<dim>(n_points_1D);
-    else if (reference_cell == Type::Wedge)
+    else if (*this == Type::Wedge)
       return Simplex::QGaussWedge<dim>(n_points_1D);
     else
       Assert(false, ExcNotImplemented());
@@ -202,10 +201,13 @@ namespace ReferenceCell
 
   template <int dim>
   const Quadrature<dim> &
-  get_nodal_type_quadrature(const Type &reference_cell)
+  Type::get_nodal_type_quadrature() const
   {
-    AssertDimension(dim, reference_cell.get_dimension());
+    AssertDimension(dim, get_dimension());
 
+    // A function that is used to fill a quadrature object of the
+    // desired type the first time we encounter a particular
+    // reference cell
     const auto create_quadrature = [](const Type &reference_cell) {
       Triangulation<dim> tria;
       GridGenerator::reference_cell(reference_cell, tria);
@@ -213,35 +215,30 @@ namespace ReferenceCell
       return Quadrature<dim>(tria.get_vertices());
     };
 
-    if (reference_cell == Type::get_hypercube<dim>())
+    if (is_hyper_cube())
       {
-        static const Quadrature<dim> quadrature =
-          create_quadrature(reference_cell);
+        static const Quadrature<dim> quadrature = create_quadrature(*this);
         return quadrature;
       }
-    else if (reference_cell == Type::Tri || reference_cell == Type::Tet)
+    else if (is_simplex())
       {
-        static const Quadrature<dim> quadrature =
-          create_quadrature(reference_cell);
+        static const Quadrature<dim> quadrature = create_quadrature(*this);
         return quadrature;
       }
-    else if (reference_cell == Type::Pyramid)
+    else if (*this == Type::Pyramid)
       {
-        static const Quadrature<dim> quadrature =
-          create_quadrature(reference_cell);
+        static const Quadrature<dim> quadrature = create_quadrature(*this);
         return quadrature;
       }
-    else if (reference_cell == Type::Wedge)
+    else if (*this == Type::Wedge)
       {
-        static const Quadrature<dim> quadrature =
-          create_quadrature(reference_cell);
+        static const Quadrature<dim> quadrature = create_quadrature(*this);
         return quadrature;
       }
     else
       Assert(false, ExcNotImplemented());
 
-    static Quadrature<dim> dummy;
-
+    static const Quadrature<dim> dummy;
     return dummy; // never reached
   }
 
index 418844ce498402348a291e8da5cc5f7c47e52b23..be85dfd4f113a15ed6d611e1400feb5fa09ff6ee 100644 (file)
@@ -32,8 +32,8 @@ for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS)
 
 for (deal_II_dimension : DIMENSIONS)
   {
-    template Quadrature<deal_II_dimension> get_gauss_type_quadrature(
-      const Type &reference_cell, const unsigned n_points_1D);
-    template const Quadrature<deal_II_dimension> &get_nodal_type_quadrature(
-      const Type &reference_cell);
+    template Quadrature<deal_II_dimension> Type::get_gauss_type_quadrature(
+      const unsigned n_points_1D) const;
+    template const Quadrature<deal_II_dimension>
+      &Type::get_nodal_type_quadrature() const;
   }
index 9a45e95ed85dca62ef032695dc5e21ac277d5ec8..3cac6f262f15593ca7c98c841c5d12ec337e0d00 100644 (file)
@@ -47,7 +47,7 @@ compute_nodal_quadrature(const FiniteElement<dim, spacedim> &fe)
   const ReferenceCell::Type type = fe.reference_cell_type();
 
   const Quadrature<dim> q_gauss =
-    ReferenceCell::get_gauss_type_quadrature<dim>(type, fe.tensor_degree() + 1);
+    type.get_gauss_type_quadrature<dim>(fe.tensor_degree() + 1);
   Triangulation<dim, spacedim> tria;
   GridGenerator::reference_cell(type, tria);
   const Mapping<dim, spacedim> &mapping =

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.