]> https://gitweb.dealii.org/ - dealii.git/commitdiff
QProjector: combine triangle and quadrilateral code. 18306/head
authorDavid Wells <drwells@email.unc.edu>
Mon, 31 Mar 2025 13:17:15 +0000 (09:17 -0400)
committerDavid Wells <drwells@email.unc.edu>
Tue, 1 Apr 2025 11:00:06 +0000 (07:00 -0400)
This also fixes some bugs where, in mixed or periodic meshes, we may
need quadrature rules defined on quadrilateral faces in the reversed
orientation (which we previously assumed did not exist).

12 files changed:
source/base/qprojector.cc
source/fe/fe_abf.cc
source/fe/fe_bdm.cc
source/fe/fe_poly_tensor.cc
source/fe/fe_raviart_thomas.cc
source/fe/fe_rt_bubbles.cc
tests/base/qprojector.output
tests/base/quadrature_sorted_test.cc
tests/base/quadrature_test.cc
tests/base/quadrature_test.output
tests/fe/rt_normal_02.output
tests/simplex/variable_face_quadratures_01.cc

index ee9ed6aac17ac632fca9461553930fc2aa7b564f..bfd32ae99b7cbdd94ddda3783fddeb98cdce70e6 100644 (file)
@@ -758,122 +758,53 @@ Quadrature<2>
 QProjector<2>::project_to_all_faces(const ReferenceCell      &reference_cell,
                                     const hp::QCollection<1> &quadrature)
 {
-  if (reference_cell == ReferenceCells::Triangle)
-    {
-      const auto support_points_line =
-        [](const auto &face, const auto &orientation) -> std::vector<Point<2>> {
-        // MSVC struggles when using face.first.begin()
-        const Point<2, double>   *vertices_ptr = &face.first[0];
-        ArrayView<const Point<2>> vertices(vertices_ptr, face.first.size());
-        const auto                temp =
-          ReferenceCells::Line.permute_by_combined_orientation(vertices,
-                                                               orientation);
-        return std::vector<Point<2>>(temp.begin(),
-                                     temp.begin() + face.first.size());
-      };
-
-      // reference faces (defined by its support points and arc length)
-      const std::array<std::pair<std::array<Point<2>, 2>, double>, 3> faces = {
-        {{{{Point<2>(0.0, 0.0), Point<2>(1.0, 0.0)}}, 1.0},
-         {{{Point<2>(1.0, 0.0), Point<2>(0.0, 1.0)}}, std::sqrt(2.0)},
-         {{{Point<2>(0.0, 1.0), Point<2>(0.0, 0.0)}}, 1.0}}};
-
-      // linear polynomial to map the reference quadrature points correctly
-      // on faces
-      const auto poly = BarycentricPolynomials<1>::get_fe_p_basis(1);
-
-      // new (projected) quadrature points and weights
-      std::vector<Point<2>> points;
-      std::vector<double>   weights;
-
-      // loop over all faces (lines) ...
-      for (unsigned int face_no = 0; face_no < faces.size(); ++face_no)
-        // ... and over all possible orientations
-        for (unsigned int orientation = 0; orientation < 2; ++orientation)
-          {
-            const auto &face = faces[face_no];
-
-            // determine support point of the current line with the correct
-            // orientation
-            std::vector<Point<2>> support_points =
-              support_points_line(face, orientation);
+  // new (projected) quadrature points and weights
+  std::vector<Point<2>> points;
+  std::vector<double>   weights;
+
+  // loop over all faces (lines) ...
+  for (const unsigned int face_no : reference_cell.face_indices())
+    // ... and over all possible orientations
+    for (types::geometric_orientation orientation = 0;
+         orientation < reference_cell.n_face_orientations(face_no);
+         ++orientation)
+      {
+        std::array<Point<2>, 2> support_points{
+          {reference_cell.face_vertex_location<2>(face_no, 0),
+           reference_cell.face_vertex_location<2>(face_no, 1)}};
+        Assert(orientation == numbers::default_geometric_orientation ||
+                 orientation == numbers::reverse_line_orientation,
+               ExcInternalError());
+        if (orientation == numbers::reverse_line_orientation)
+          std::swap(support_points[0], support_points[1]);
 
-            // the quadrature rule to be projected ...
-            const auto &sub_quadrature_points =
-              quadrature[quadrature.size() == 1 ? 0 : face_no].get_points();
-            const auto &sub_quadrature_weights =
-              quadrature[quadrature.size() == 1 ? 0 : face_no].get_weights();
+        // the quadrature rule to be projected ...
+        const auto &sub_quadrature_points =
+          quadrature[quadrature.size() == 1 ? 0 : face_no].get_points();
+        const auto &sub_quadrature_weights =
+          quadrature[quadrature.size() == 1 ? 0 : face_no].get_weights();
 
-            // loop over all quadrature points
-            for (unsigned int j = 0; j < sub_quadrature_points.size(); ++j)
-              {
-                Point<2> mapped_point;
+        // loop over all quadrature points
+        for (unsigned int j = 0; j < sub_quadrature_points.size(); ++j)
+          {
+            Point<2> mapped_point;
 
-                // map reference quadrature point
-                for (unsigned int i = 0; i < 2; ++i)
-                  mapped_point +=
-                    support_points[i] *
-                    poly.compute_value(i, sub_quadrature_points[j]);
+            // map reference quadrature point
+            for (unsigned int i = 0; i < support_points.size(); ++i)
+              mapped_point +=
+                support_points[i] *
+                ReferenceCells::Line.template d_linear_shape_function<1>(
+                  sub_quadrature_points[j], i);
 
-                points.emplace_back(mapped_point);
+            points.emplace_back(mapped_point);
 
-                // scale weight by arc length
-                weights.emplace_back(sub_quadrature_weights[j] * face.second);
-              }
+            // scale weight by arc length
+            weights.emplace_back(sub_quadrature_weights[j] *
+                                 reference_cell.face_measure(face_no));
           }
+      }
 
-      // construct new quadrature rule
-      return Quadrature<2>(std::move(points), std::move(weights));
-    }
-
-  Assert(reference_cell == ReferenceCells::Quadrilateral, ExcNotImplemented());
-
-  const unsigned int dim = 2;
-
-  const unsigned int n_faces = GeometryInfo<dim>::faces_per_cell;
-
-  unsigned int n_points_total = 0;
-
-  if (quadrature.size() == 1)
-    n_points_total = quadrature[0].size() * GeometryInfo<dim>::faces_per_cell;
-  else
-    {
-      AssertDimension(quadrature.size(), GeometryInfo<dim>::faces_per_cell);
-      for (const auto &q : quadrature)
-        n_points_total += q.size();
-    }
-
-  // first fix quadrature points
-  std::vector<Point<dim>> q_points;
-  q_points.reserve(n_points_total);
-  std::vector<Point<dim>> help;
-  help.reserve(quadrature.max_n_quadrature_points());
-
-  // project to each face and append
-  // results
-  for (unsigned int face = 0; face < n_faces; ++face)
-    {
-      help.resize(quadrature[quadrature.size() == 1 ? 0 : face].size());
-      project_to_face(reference_cell,
-                      quadrature[quadrature.size() == 1 ? 0 : face],
-                      face,
-                      help);
-      std::copy(help.begin(), help.end(), std::back_inserter(q_points));
-    }
-
-  // next copy over weights
-  std::vector<double> weights;
-  weights.reserve(n_points_total);
-  for (unsigned int face = 0; face < n_faces; ++face)
-    std::copy(
-      quadrature[quadrature.size() == 1 ? 0 : face].get_weights().begin(),
-      quadrature[quadrature.size() == 1 ? 0 : face].get_weights().end(),
-      std::back_inserter(weights));
-
-  Assert(q_points.size() == n_points_total, ExcInternalError());
-  Assert(weights.size() == n_points_total, ExcInternalError());
-
-  return Quadrature<dim>(std::move(q_points), std::move(weights));
+  return Quadrature<2>(std::move(points), std::move(weights));
 }
 
 
@@ -1361,41 +1292,20 @@ QProjector<dim>::DataSetDescriptor::face(
   Assert(dim == 1 ||
            (combined_orientation < reference_cell.n_face_orientations(face_no)),
          ExcInternalError());
-  if (reference_cell == ReferenceCells::Triangle ||
-      reference_cell == ReferenceCells::Tetrahedron)
-    {
-      if (dim == 2)
-        return {(2 * face_no + (combined_orientation ==
-                                    numbers::default_geometric_orientation ?
-                                  1 :
-                                  0)) *
-                n_quadrature_points};
-      else if (dim == 3)
-        {
-          return {(reference_cell.n_face_orientations(face_no) * face_no +
-                   combined_orientation) *
-                  n_quadrature_points};
-        }
-    }
 
-  Assert(reference_cell == ReferenceCells::get_hypercube<dim>(),
-         ExcNotImplemented());
-
-  Assert(face_no < GeometryInfo<dim>::faces_per_cell, ExcInternalError());
-
-  switch (dim)
-    {
-      case 1:
-      case 2:
-        return face_no * n_quadrature_points;
-      case 3:
-        return (face_no +
-                GeometryInfo<dim>::faces_per_cell * combined_orientation) *
-               n_quadrature_points;
-      default:
-        DEAL_II_ASSERT_UNREACHABLE();
-    }
-  return numbers::invalid_unsigned_int;
+  // TODO: once the default orientation is 0 we can combine this with the
+  // general branch
+  if (reference_cell == ReferenceCells::Line)
+    return {face_no};
+  // TODO: index Hexahedra in the same way as everything else
+  else if (reference_cell == ReferenceCells::Hexahedron)
+    return (face_no +
+            GeometryInfo<dim>::faces_per_cell * combined_orientation) *
+           n_quadrature_points;
+  else
+    return {(reference_cell.n_face_orientations(face_no) * face_no +
+             combined_orientation) *
+            n_quadrature_points};
 }
 
 
@@ -1428,6 +1338,12 @@ QProjector<dim>::DataSetDescriptor::face(
   const types::geometric_orientation combined_orientation,
   const hp::QCollection<dim - 1>    &quadrature)
 {
+  // TODO: once we move to representing the default orientation as 0 (instead of
+  // 1) we can get rid of the dim = 1 check
+  Assert(dim == 1 ||
+           (combined_orientation < reference_cell.n_face_orientations(face_no)),
+         ExcInternalError());
+
   if (reference_cell == ReferenceCells::Triangle ||
       reference_cell == ReferenceCells::Tetrahedron ||
       reference_cell == ReferenceCells::Wedge ||
@@ -1459,10 +1375,9 @@ QProjector<dim>::DataSetDescriptor::face(
           offset += scale[i] * quadrature[i].size();
 
       if (dim == 2)
-        return {
-          offset +
-          (combined_orientation == numbers::default_geometric_orientation) *
-            quadrature[quadrature.size() == 1 ? 0 : face_no].size()};
+        return {offset +
+                combined_orientation *
+                  quadrature[quadrature.size() == 1 ? 0 : face_no].size()};
       else if (dim == 3)
         {
           return {offset +
@@ -1479,17 +1394,22 @@ QProjector<dim>::DataSetDescriptor::face(
   switch (dim)
     {
       case 1:
+        return face_no;
       case 2:
         {
+          unsigned int offset = 0;
+
           if (quadrature.size() == 1)
-            return quadrature[0].size() * face_no;
+            offset = reference_cell.n_face_orientations(0) *
+                     quadrature[0].size() * face_no;
           else
-            {
-              unsigned int result = 0;
-              for (unsigned int i = 0; i < face_no; ++i)
-                result += quadrature[i].size();
-              return result;
-            }
+            for (unsigned int i = 0; i < face_no; ++i)
+              offset += reference_cell.n_face_orientations(face_no) *
+                        quadrature[i].size();
+
+          return {offset +
+                  combined_orientation *
+                    quadrature[quadrature.size() == 1 ? 0 : face_no].size()};
         }
       case 3:
         {
index be10e3f6fc947b54c8c951bce4b6a41a9e1dc3f6..98bce41a77897d40b860ef22891cbb5a7445c2d6 100644 (file)
@@ -234,12 +234,23 @@ FE_ABF<dim>::initialize_support_points(const unsigned int deg)
       Quadrature<dim> faces =
         QProjector<dim>::project_to_all_faces(this->reference_cell(),
                                               face_points);
-      for (; current < GeometryInfo<dim>::faces_per_cell * n_face_points;
-           ++current)
+      for (unsigned int face_no = 0;
+           face_no < GeometryInfo<dim>::faces_per_cell;
+           ++face_no)
         {
-          // Enter the support point
-          // into the vector
-          this->generalized_support_points[current] = faces.point(current);
+          const auto offset = QProjector<dim>::DataSetDescriptor::face(
+            this->reference_cell(),
+            face_no,
+            numbers::default_geometric_orientation,
+            n_face_points);
+          for (unsigned int face_point = 0; face_point < n_face_points;
+               ++face_point)
+            {
+              // Enter the support point into the vector
+              this->generalized_support_points[current] =
+                faces.point(offset + face_point);
+              ++current;
+            }
         }
 
 
index f52a220f4e9ea4e16430802ee56d9d79774d3ca3..e89d8693d4ac60d8abe7a5f13a3b2b74bbd3b100 100644 (file)
@@ -387,16 +387,19 @@ FE_BDM<dim>::initialize_support_points(const unsigned int deg)
   const Quadrature<dim> faces =
     QProjector<dim>::project_to_all_faces(this->reference_cell(), face_points);
 
-  for (unsigned int f = 0; f < GeometryInfo<dim>::faces_per_cell; ++f)
+  for (unsigned int face_no = 0; face_no < GeometryInfo<dim>::faces_per_cell;
+       ++face_no)
     {
       const auto offset = QProjector<dim>::DataSetDescriptor::face(
         this->reference_cell(),
-        f,
+        face_no,
         numbers::default_geometric_orientation,
         face_points.size());
-      for (unsigned int k = 0; k < face_points.size(); ++k)
-        this->generalized_support_points[face_points.size() * f + k] =
-          faces.point(offset + k);
+      for (unsigned int face_point = 0; face_point < face_points.size();
+           ++face_point)
+        this->generalized_support_points[face_points.size() * face_no +
+                                         face_point] =
+          faces.point(offset + face_point);
     }
 
   // Currently, for backward compatibility, we do not use moments, but
index f17a13ee80e7c043e32cd6ac7b7aeab1395d8887..3efe77137d5a809536e831ab8399d7a847adc2c2 100644 (file)
@@ -1180,12 +1180,21 @@ FE_PolyTensor<dim, spacedim>::fill_fe_face_values(
   // to take (all data sets for all
   // faces are stored contiguously)
 
-  const auto offset =
-    QProjector<dim>::DataSetDescriptor::face(this->reference_cell(),
-                                             face_no,
-                                             cell->combined_face_orientation(
-                                               face_no),
-                                             n_q_points);
+  // TODO: The same 'legacy' comments for 2d apply here as well: these classes
+  // do not handle non-standard orientations in 2d in a way consistent with the
+  // rest of the library, but are consistent with themselves (see, e.g., the
+  // fe_conformity_dim_2 tests).
+  //
+  // In this case: all of this code was written assuming that QProjector assumed
+  // that all faces were in the default orientation in 2d, but contains special
+  // workarounds in case that isn't the case. Hence, to keep those workarounds
+  // working, we still assume that all faces are in the default orientation.
+  const auto offset = QProjector<dim>::DataSetDescriptor::face(
+    this->reference_cell(),
+    face_no,
+    dim == 2 ? numbers::default_geometric_orientation :
+               cell->combined_face_orientation(face_no),
+    n_q_points);
 
   // TODO: Size assertions
 
index 378e984d52bc8ce6a2d6e7f0a446d6cc83fc0745..40719719665878cf165b2220e646b98a8463c8e9 100644 (file)
@@ -208,16 +208,23 @@ FE_RaviartThomas<dim>::initialize_support_points(const unsigned int deg)
         QProjector<dim>::project_to_all_faces(this->reference_cell(),
                                               face_points);
 
-      for (; current < this->reference_cell().n_faces() * n_face_points;
-           ++current)
+      for (unsigned int face_no = 0;
+           face_no < GeometryInfo<dim>::faces_per_cell;
+           ++face_no)
         {
-          // Enter the support point into the vector
-          this->generalized_support_points[current] =
-            faces.point(current + QProjector<dim>::DataSetDescriptor::face(
-                                    this->reference_cell(),
-                                    0,
-                                    numbers::default_geometric_orientation,
-                                    n_face_points));
+          const auto offset = QProjector<dim>::DataSetDescriptor::face(
+            this->reference_cell(),
+            face_no,
+            numbers::default_geometric_orientation,
+            n_face_points);
+          for (unsigned int face_point = 0; face_point < n_face_points;
+               ++face_point)
+            {
+              // Enter the support point into the vector
+              this->generalized_support_points[current] =
+                faces.point(offset + face_point);
+              ++current;
+            }
         }
     }
 
index 1f9769b8c67ec8ec13f67d28d56bcbc1d42f9375..7349e54d508af2fcd3d29be5787650d2ac9303dd 100644 (file)
@@ -184,18 +184,24 @@ FE_RT_Bubbles<dim>::initialize_support_points(const unsigned int deg)
       Quadrature<dim> faces =
         QProjector<dim>::project_to_all_faces(this->reference_cell(),
                                               face_points);
-      for (unsigned int k = 0; k < this->n_dofs_per_face(face_no) *
-                                     GeometryInfo<dim>::faces_per_cell;
-           ++k)
-        this->generalized_support_points[k] =
-          faces.point(k + QProjector<dim>::DataSetDescriptor::face(
-                            this->reference_cell(),
-                            0,
-                            numbers::default_geometric_orientation,
-                            this->n_dofs_per_face(face_no)));
-
-      current =
-        this->n_dofs_per_face(face_no) * GeometryInfo<dim>::faces_per_cell;
+      for (unsigned int face_no = 0;
+           face_no < GeometryInfo<dim>::faces_per_cell;
+           ++face_no)
+        {
+          const auto offset = QProjector<dim>::DataSetDescriptor::face(
+            this->reference_cell(),
+            face_no,
+            numbers::default_geometric_orientation,
+            face_points.size());
+          for (unsigned int face_point = 0; face_point < face_points.size();
+               ++face_point)
+            {
+              // Enter the support point into the vector
+              this->generalized_support_points[current] =
+                faces.point(offset + face_point);
+              ++current;
+            }
+        }
     }
 
   if (deg == 1)
index 94997d10cf3e5803d958ff5d5dbfe8648a9c23f3..4b8265450265f2b620da67fa798acf00f73c716c 100644 (file)
@@ -1,8 +1,8 @@
 
 DEAL::
-DEAL:line::length: 0
-DEAL:line::length: 0
-DEAL:line::length: 0
+DEAL:line::length: 0.0
+DEAL:line::length: 0.0
+DEAL:line::length: 0.0
 DEAL:face::Checking dim 1 1d-points 0
 DEAL:face::Face 0
 DEAL:face::0.0
@@ -431,26 +431,26 @@ DEAL:face::0.0 1.0 1.0
 DEAL:face::1.0 1.0 1.0
 DEAL:all::Checking dim 2 1d-points 2
 DEAL:all::Face 0 orientation false
-DEAL:all::0.0 0.0
 DEAL:all::0.0 1.0
+DEAL:all::0.0 0.0
 DEAL:all::Face 0 orientation true
 DEAL:all::0.0 0.0
 DEAL:all::0.0 1.0
 DEAL:all::Face 1 orientation false
-DEAL:all::1.0 0.0
 DEAL:all::1.0 1.0
+DEAL:all::1.0 0.0
 DEAL:all::Face 1 orientation true
 DEAL:all::1.0 0.0
 DEAL:all::1.0 1.0
 DEAL:all::Face 2 orientation false
-DEAL:all::0.0 0.0
 DEAL:all::1.0 0.0
+DEAL:all::0.0 0.0
 DEAL:all::Face 2 orientation true
 DEAL:all::0.0 0.0
 DEAL:all::1.0 0.0
 DEAL:all::Face 3 orientation false
-DEAL:all::0.0 1.0
 DEAL:all::1.0 1.0
+DEAL:all::0.0 1.0
 DEAL:all::Face 3 orientation true
 DEAL:all::0.0 1.0
 DEAL:all::1.0 1.0
@@ -889,33 +889,33 @@ DEAL:face::0.50 1.0 1.0
 DEAL:face::1.0 1.0 1.0
 DEAL:all::Checking dim 2 1d-points 3
 DEAL:all::Face 0 orientation false
-DEAL:all::0.0 0.0
-DEAL:all::0.0 0.50
 DEAL:all::0.0 1.0
+DEAL:all::0.0 0.50
+DEAL:all::0.0 0.0
 DEAL:all::Face 0 orientation true
 DEAL:all::0.0 0.0
 DEAL:all::0.0 0.50
 DEAL:all::0.0 1.0
 DEAL:all::Face 1 orientation false
-DEAL:all::1.0 0.0
-DEAL:all::1.0 0.50
 DEAL:all::1.0 1.0
+DEAL:all::1.0 0.50
+DEAL:all::1.0 0.0
 DEAL:all::Face 1 orientation true
 DEAL:all::1.0 0.0
 DEAL:all::1.0 0.50
 DEAL:all::1.0 1.0
 DEAL:all::Face 2 orientation false
-DEAL:all::0.0 0.0
-DEAL:all::0.50 0.0
 DEAL:all::1.0 0.0
+DEAL:all::0.50 0.0
+DEAL:all::0.0 0.0
 DEAL:all::Face 2 orientation true
 DEAL:all::0.0 0.0
 DEAL:all::0.50 0.0
 DEAL:all::1.0 0.0
 DEAL:all::Face 3 orientation false
-DEAL:all::0.0 1.0
-DEAL:all::0.50 1.0
 DEAL:all::1.0 1.0
+DEAL:all::0.50 1.0
+DEAL:all::0.0 1.0
 DEAL:all::Face 3 orientation true
 DEAL:all::0.0 1.0
 DEAL:all::0.50 1.0
@@ -1925,11 +1925,11 @@ DEAL:face::0.75 1.0 1.0
 DEAL:face::1.0 1.0 1.0
 DEAL:all::Checking dim 2 1d-points 5
 DEAL:all::Face 0 orientation false
-DEAL:all::0.0 0.0
-DEAL:all::0.0 0.25
-DEAL:all::0.0 0.50
-DEAL:all::0.0 0.75
 DEAL:all::0.0 1.0
+DEAL:all::0.0 0.75
+DEAL:all::0.0 0.50
+DEAL:all::0.0 0.25
+DEAL:all::0.0 0.0
 DEAL:all::Face 0 orientation true
 DEAL:all::0.0 0.0
 DEAL:all::0.0 0.25
@@ -1937,11 +1937,11 @@ DEAL:all::0.0 0.50
 DEAL:all::0.0 0.75
 DEAL:all::0.0 1.0
 DEAL:all::Face 1 orientation false
-DEAL:all::1.0 0.0
-DEAL:all::1.0 0.25
-DEAL:all::1.0 0.50
-DEAL:all::1.0 0.75
 DEAL:all::1.0 1.0
+DEAL:all::1.0 0.75
+DEAL:all::1.0 0.50
+DEAL:all::1.0 0.25
+DEAL:all::1.0 0.0
 DEAL:all::Face 1 orientation true
 DEAL:all::1.0 0.0
 DEAL:all::1.0 0.25
@@ -1949,11 +1949,11 @@ DEAL:all::1.0 0.50
 DEAL:all::1.0 0.75
 DEAL:all::1.0 1.0
 DEAL:all::Face 2 orientation false
-DEAL:all::0.0 0.0
-DEAL:all::0.25 0.0
-DEAL:all::0.50 0.0
-DEAL:all::0.75 0.0
 DEAL:all::1.0 0.0
+DEAL:all::0.75 0.0
+DEAL:all::0.50 0.0
+DEAL:all::0.25 0.0
+DEAL:all::0.0 0.0
 DEAL:all::Face 2 orientation true
 DEAL:all::0.0 0.0
 DEAL:all::0.25 0.0
@@ -1961,11 +1961,11 @@ DEAL:all::0.50 0.0
 DEAL:all::0.75 0.0
 DEAL:all::1.0 0.0
 DEAL:all::Face 3 orientation false
-DEAL:all::0.0 1.0
-DEAL:all::0.25 1.0
-DEAL:all::0.50 1.0
-DEAL:all::0.75 1.0
 DEAL:all::1.0 1.0
+DEAL:all::0.75 1.0
+DEAL:all::0.50 1.0
+DEAL:all::0.25 1.0
+DEAL:all::0.0 1.0
 DEAL:all::Face 3 orientation true
 DEAL:all::0.0 1.0
 DEAL:all::0.25 1.0
index dfaaa3151925f00519aa95e160f0c40a1ecf1f69..3257d14fb20fb5cf114006640c3ca6b4a9d7d035 100644 (file)
@@ -155,7 +155,11 @@ check_faces(const std::vector<Quadrature<dim - 1> *> &quadratures,
           switch (dim)
             {
               case 2:
-                exact_int = 2 * (sub ? 2 : 1) / (double)(i + 1);
+                // TODO: once we support multiple orientations per subface in 2d
+                // we can get rid of the second ternary if here (because we will
+                // always integrate over all orientations: i.e., both of them in
+                // 2d).
+                exact_int = 2 * (sub ? 2 : 1) * (sub ? 1 : 2) / (double)(i + 1);
                 break;
               case 3:
                 exact_int =
@@ -170,7 +174,10 @@ check_faces(const std::vector<Quadrature<dim - 1> *> &quadratures,
       // over the whole surface (all
       // combinations of face_orientation,
       // face_flip and face_rotation)
-      while (err < (dim == 3 ? 8 : 1) * 2e-14);
+      // but 2 in 2d (only two possible
+      // orientations) when we are not using
+      // subfaces (same as the previous comment)
+      while (err < (dim == 3 ? 8 : (sub ? 1 : 2)) * 2e-14);
       // Uncomment here for testing
       //      deallog << " (Int " << quadrature_int << '-' << exact_int << '='
       //      << err << ')';
index d51e64f1d5e191542cce5a83feb9fea35b79884e..6009fb82e632bc16e03185a839ffde7fe646ca9d 100644 (file)
@@ -176,7 +176,11 @@ check_faces(const std::vector<Quadrature<dim - 1> *> &quadratures,
           switch (dim)
             {
               case 2:
-                exact_int = 2 * (sub ? 2 : 1) / (double)(i + 1);
+                // TODO: once we support multiple orientations per subface in 2d
+                // we can get rid of the second ternary if here (because we will
+                // always integrate over all orientations: i.e., both of them in
+                // 2d).
+                exact_int = 2 * (sub ? 2 : 1) * (sub ? 1 : 2) / (double)(i + 1);
                 break;
               case 3:
                 exact_int =
@@ -191,7 +195,10 @@ check_faces(const std::vector<Quadrature<dim - 1> *> &quadratures,
       // over the whole surface (all
       // combinations of face_orientation,
       // face_flip and face_rotation)
-      while (err < (dim == 3 ? 8 : 1) * 2e-14);
+      // but 2 in 2d (only two possible
+      // orientations) when we are not using
+      // subfaces (same as the previous comment)
+      while (err < (dim == 3 ? 8 : (sub ? 1 : 2)) * 2e-14);
       // Uncomment here for testing
       //      deallog << " (Int " << quadrature_int << '-' << exact_int << '='
       //      << err << ')';
index 42383d7cf9b3824cec8ec7871e8fd228302811e6..1dc3049d283c3a704f0ceed27856cfa0c1693c0f 100644 (file)
@@ -94,22 +94,22 @@ DEAL:2d:faces::Quadrature no.17 is exact for polynomials of degree 5
 DEAL:2d:faces::Quadrature no.18 is exact for polynomials of degree 7
 DEAL:2d:faces::Quadrature no.19 is exact for polynomials of degree 9
 DEAL:2d:faces::Quadrature no.20 is exact for polynomials of degree 11
-DEAL:2d:faces::Quadrature no.21 is exact for polynomials of degree 0
-DEAL:2d:faces::Quadrature no.22 is exact for polynomials of degree 0
-DEAL:2d:faces::Quadrature no.23 is exact for polynomials of degree 2
-DEAL:2d:faces::Quadrature no.24 is exact for polynomials of degree 2
-DEAL:2d:faces::Quadrature no.25 is exact for polynomials of degree 4
-DEAL:2d:faces::Quadrature no.26 is exact for polynomials of degree 4
-DEAL:2d:faces::Quadrature no.27 is exact for polynomials of degree 6
-DEAL:2d:faces::Quadrature no.28 is exact for polynomials of degree 6
-DEAL:2d:faces::Quadrature no.29 is exact for polynomials of degree 8
-DEAL:2d:faces::Quadrature no.30 is exact for polynomials of degree 8
-DEAL:2d:faces::Quadrature no.31 is exact for polynomials of degree 10
-DEAL:2d:faces::Quadrature no.32 is exact for polynomials of degree 10
-DEAL:2d:faces::Quadrature no.33 is exact for polynomials of degree 12
-DEAL:2d:faces::Quadrature no.34 is exact for polynomials of degree 12
-DEAL:2d:faces::Quadrature no.35 is exact for polynomials of degree 14
-DEAL:2d:faces::Quadrature no.36 is exact for polynomials of degree 14
+DEAL:2d:faces::Quadrature no.21 is exact for polynomials of degree 1
+DEAL:2d:faces::Quadrature no.22 is exact for polynomials of degree 1
+DEAL:2d:faces::Quadrature no.23 is exact for polynomials of degree 3
+DEAL:2d:faces::Quadrature no.24 is exact for polynomials of degree 3
+DEAL:2d:faces::Quadrature no.25 is exact for polynomials of degree 5
+DEAL:2d:faces::Quadrature no.26 is exact for polynomials of degree 5
+DEAL:2d:faces::Quadrature no.27 is exact for polynomials of degree 7
+DEAL:2d:faces::Quadrature no.28 is exact for polynomials of degree 7
+DEAL:2d:faces::Quadrature no.29 is exact for polynomials of degree 9
+DEAL:2d:faces::Quadrature no.30 is exact for polynomials of degree 9
+DEAL:2d:faces::Quadrature no.31 is exact for polynomials of degree 11
+DEAL:2d:faces::Quadrature no.32 is exact for polynomials of degree 11
+DEAL:2d:faces::Quadrature no.33 is exact for polynomials of degree 13
+DEAL:2d:faces::Quadrature no.34 is exact for polynomials of degree 13
+DEAL:2d:faces::Quadrature no.35 is exact for polynomials of degree 15
+DEAL:2d:faces::Quadrature no.36 is exact for polynomials of degree 15
 DEAL:2d:subfaces::Quadrature no.0 is exact for polynomials of degree 1
 DEAL:2d:subfaces::Quadrature no.1 is exact for polynomials of degree 1
 DEAL:2d:subfaces::Quadrature no.2 is exact for polynomials of degree 3
index 8ae862c75d4269467f62c1d589be9b5cc966c570..655070f0410282cd9842edd3bf457dbcf0decd0a 100644 (file)
@@ -1,5 +1,5 @@
 
-DEAL::n_q_face=6, n_q_proj=24
+DEAL::n_q_face=6, n_q_proj=48
 DEAL::Testing face with center at 0.24866160 0.03044371
 DEAL::  QP=0, error=0.00000000, u.n=-8.30976642, u_neighbor.n=-8.30976642
 DEAL::  QP=1, error=0.00000000, u.n=-2.15456843, u_neighbor.n=-2.15456843
index 3005c14207193dc8f10570a60264dd84c4ca765c..3d972b927ad3cea98c41184decc93a7c8c4b94b1 100644 (file)
@@ -83,9 +83,7 @@ test<2>()
                         i = QProjector<dim>::DataSetDescriptor::face(
                           ReferenceCells::Quadrilateral,
                           face_no,
-                          false,
-                          false,
-                          false,
+                          numbers::default_geometric_orientation,
                           quad_ref);
            q < quad_ref[face_no].size();
            ++q, ++i)

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.