]> https://gitweb.dealii.org/ - dealii.git/commitdiff
QProjector: consistently call the orientation the combined_orientation. 18374/head
authorDavid Wells <drwells@email.unc.edu>
Mon, 21 Apr 2025 13:24:44 +0000 (09:24 -0400)
committerDavid Wells <drwells@email.unc.edu>
Mon, 21 Apr 2025 20:34:25 +0000 (16:34 -0400)
include/deal.II/base/qprojector.h
source/base/qprojector.cc

index 269d744d2ea00ec9f8f308600e4b84f7f7f08b3d..91a26b12679ab679a4d08a2c620f56bd6a1594ab 100644 (file)
@@ -58,11 +58,11 @@ class ReferenceCell;
  *
  * In practice, computing face integrals (e.g., via FEFaceValues or
  * FESubfaceValues) requires quadrature rules for all possible permutations of
- * face number, subface number, and orientation. This class provides several
- * functions for doing just that, such as QProjector::project_to_all_faces().
- * Furthermore, the DataSetDescriptor class implements indexing for converting
- * between face number, subface number, and orientation to the index of the
- * associated Quadrature rule.
+ * face number, subface number, and combined orientation. This class provides
+ * several functions for doing just that, such as
+ * QProjector::project_to_all_faces(). Furthermore, the DataSetDescriptor class
+ * implements indexing for converting between face number, subface number, and
+ * combined orientation to the index of the associated Quadrature rule.
  */
 template <int dim>
 class QProjector
@@ -136,7 +136,7 @@ public:
   project_to_face(const ReferenceCell               &reference_cell,
                   const SubQuadrature               &quadrature,
                   const unsigned int                 face_no,
-                  const types::geometric_orientation orientation);
+                  const types::geometric_orientation combined_orientation);
 
   /**
    * Compute the quadrature points on the cell if the given quadrature formula
@@ -147,12 +147,12 @@ public:
    * @note Only the points are transformed. The quadrature weights are the
    * same as those of the original rule.
    *
-   * @deprecated Use the version of project_to_subface() which takes an
-   * orientation argument instead.
+   * @deprecated Use the version of project_to_subface() which takes a
+   * combined_orientation argument instead.
    */
   DEAL_II_DEPRECATED_EARLY_WITH_COMMENT(
-    "Use the version of project_to_subface() which takes an orientation "
-    "argument instead.")
+    "Use the version of project_to_subface() which takes a "
+    "combined_orientation argument instead.")
   static void
   project_to_subface(const ReferenceCell           &reference_cell,
                      const SubQuadrature           &quadrature,
@@ -175,12 +175,12 @@ public:
    * that the cell is a line (1D), a quad (2d), or a hex (3d). Use the other
    * version of this function that takes the reference cell type instead.
    *
-   * @deprecated Use the version of project_to_subface() which takes an
-   * orientation argument instead.
+   * @deprecated Use the version of project_to_subface() which takes a
+   * combined_orientation argument instead.
    */
   DEAL_II_DEPRECATED_EARLY_WITH_COMMENT(
-    "Use the version of project_to_subface() which takes an orientation "
-    "argument instead.")
+    "Use the version of project_to_subface() which takes a "
+    "combined_orientation argument instead.")
   static Quadrature<dim>
   project_to_subface(const ReferenceCell           &reference_cell,
                      const SubQuadrature           &quadrature,
@@ -198,12 +198,12 @@ public:
    * @note Only the points are transformed. The quadrature weights are the
    * same as those of the original rule.
    *
-   * @deprecated Use the version of project_to_subface() which takes an
-   * orientation argument instead.
+   * @deprecated Use the version of project_to_subface() which takes a
+   * combined_orientation argument instead.
    */
   DEAL_II_DEPRECATED_EARLY_WITH_COMMENT(
-    "Use the version of project_to_subface() which takes an orientation "
-    "argument instead.")
+    "Use the version of project_to_subface() which takes a "
+    "combined_orientation argument instead.")
   static Quadrature<dim>
   project_to_oriented_subface(const ReferenceCell             &reference_cell,
                               const SubQuadrature             &quadrature,
@@ -228,7 +228,7 @@ public:
                      const SubQuadrature               &quadrature,
                      const unsigned int                 face_no,
                      const unsigned int                 subface_no,
-                     const types::geometric_orientation orientation,
+                     const types::geometric_orientation combined_orientation,
                      const RefinementCase<dim - 1>     &ref_case);
 
   /**
@@ -245,8 +245,8 @@ public:
    * FEFaceValues.
    *
    * @note This function creates ReferenceCell::n_face_orientations() sets of
-   * quadrature points for each face which are indexed (by orientation and face
-   * number) by a DataSetDescriptor.
+   * quadrature points for each face which are indexed (by combined orientation
+   * and face number) by a DataSetDescriptor.
    */
   static Quadrature<dim>
   project_to_all_faces(const ReferenceCell            &reference_cell,
@@ -325,10 +325,10 @@ public:
    *
    * The functions QProjector::project_to_all_faces() and
    * QProjector::project_to_all_subfaces() each combine all quadrature rules
-   * (i.e., all possible combinations of face, subface, and orientation) into a
-   * single Quadrature object. DataSetDescriptor implements the correct indexing
-   * for extracting from that Quadrature rule the correct index for those
-   * values.
+   * (i.e., all possible combinations of face, subface, and combined
+   * orientation) into a single Quadrature object. DataSetDescriptor implements
+   * the correct indexing for extracting from that Quadrature rule the correct
+   * index for those values.
    */
   class DataSetDescriptor
   {
index b950161365aa8b30b26569c70d3bc4ded75682ac..2db5a4fcd914abc12709a15f58ccaee63efc8639 100644 (file)
@@ -476,13 +476,15 @@ QProjector<dim>::project_to_oriented_face(const ReferenceCell &reference_cell,
 
 template <int dim>
 Quadrature<dim>
-QProjector<dim>::project_to_face(const ReferenceCell       &reference_cell,
-                                 const Quadrature<dim - 1> &quadrature,
-                                 const unsigned int         face_no,
-                                 const types::geometric_orientation orientation)
+QProjector<dim>::project_to_face(
+  const ReferenceCell               &reference_cell,
+  const Quadrature<dim - 1>         &quadrature,
+  const unsigned int                 face_no,
+  const types::geometric_orientation combined_orientation)
 {
   AssertIndexRange(face_no, reference_cell.n_faces());
-  AssertIndexRange(orientation, reference_cell.n_face_orientations(face_no));
+  AssertIndexRange(combined_orientation,
+                   reference_cell.n_face_orientations(face_no));
   AssertDimension(reference_cell.get_dimension(), dim);
 
   std::vector<Point<dim>> points;
@@ -499,7 +501,7 @@ QProjector<dim>::project_to_face(const ReferenceCell       &reference_cell,
                                               face_vertices,
                                               reference_cell.face_measure(
                                                 face_no),
-                                              orientation,
+                                              combined_orientation,
                                               points,
                                               weights);
 
@@ -604,11 +606,12 @@ QProjector<dim>::project_to_subface(
   const SubQuadrature               &quadrature,
   const unsigned int                 face_no,
   const unsigned int                 subface_no,
-  const types::geometric_orientation orientation,
+  const types::geometric_orientation combined_orientation,
   const RefinementCase<dim - 1>     &ref_case)
 {
   AssertIndexRange(face_no, reference_cell.n_faces());
-  AssertIndexRange(orientation, reference_cell.n_face_orientations(face_no));
+  AssertIndexRange(combined_orientation,
+                   reference_cell.n_face_orientations(face_no));
   AssertDimension(reference_cell.get_dimension(), dim);
   AssertIndexRange(subface_no,
                    reference_cell.face_reference_cell(face_no)
@@ -693,7 +696,7 @@ QProjector<dim>::project_to_subface(
       else
         DEAL_II_ASSERT_UNREACHABLE();
 
-      if (orientation == numbers::reverse_line_orientation)
+      if (combined_orientation == numbers::reverse_line_orientation)
         {
           std::reverse(q_points.begin(), q_points.end());
           std::reverse(q_weights.begin(), q_weights.end());

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.