]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Use std::optional in orthogonal_equality().
authorDavid Wells <drwells@email.unc.edu>
Sun, 22 Oct 2023 01:12:03 +0000 (21:12 -0400)
committerDavid Wells <drwells@email.unc.edu>
Mon, 30 Oct 2023 20:02:19 +0000 (16:02 -0400)
This is cleaner than returning a boolean and modifying an input argument.

include/deal.II/grid/grid_tools.h
source/grid/grid_tools_dof_handlers.cc
source/grid/grid_tools_dof_handlers.inst.in
tests/dofs/dof_tools_21_b.cc
tests/dofs/dof_tools_21_c.cc

index 75c988f7a50369d4b3c17e91ddfdec171891111a..8f8eb1c09c5d1f013743f4d599130db8e542c1bd 100644 (file)
@@ -2361,6 +2361,8 @@ namespace GridTools
    *
    * @p face1 and @p face2 are considered equal, if a one to one matching
    * between its vertices can be achieved via an orthogonal equality relation.
+   * If no such relation exists then the returned std::optional object is empty
+   * (i.e., has_value() will return `false`).
    *
    * Here, two vertices <tt>v_1</tt> and <tt>v_2</tt> are considered equal, if
    * $M\cdot v_1 + offset - v_2$ is parallel to the unit vector in unit
@@ -2368,12 +2370,13 @@ namespace GridTools
    * spacedim x spacedim matrix, $M$ is set to @p matrix, otherwise $M$ is the
    * identity matrix.
    *
-   * If the matching was successful, the _relative_ orientation of @p face1
-   * with respect to @p face2 is returned in the bitset @p orientation, where
+   * If the matching was successful, the _relative_ orientation of @p face1 with
+   * respect to @p face2 is returned a std::optional<std::bitset<3>> object
+   * orientation in which
    * @code
-   * orientation[0] -> face_orientation
-   * orientation[1] -> face_flip
-   * orientation[2] -> face_rotation
+   * orientation.value()[0] = face_orientation
+   * orientation.value()[1] = face_flip
+   * orientation.value()[2] = face_rotation
    * @endcode
    *
    * In 2d <tt>face_orientation</tt> is always <tt>true</tt>,
@@ -2420,9 +2423,8 @@ namespace GridTools
    * article.
    */
   template <typename FaceIterator>
-  bool
+  std::optional<std::bitset<3>>
   orthogonal_equality(
-    std::bitset<3>                                               &orientation,
     const FaceIterator                                           &face1,
     const FaceIterator                                           &face2,
     const unsigned int                                            direction,
index 0617720f77614c0c08b4c6beae65e5498b62912c..493d088e291e340e197a3f1088aafd0c61fdc971 100644 (file)
@@ -43,6 +43,7 @@
 #include <list>
 #include <map>
 #include <numeric>
+#include <optional>
 #include <set>
 #include <vector>
 
@@ -2142,7 +2143,6 @@ namespace GridTools
     unsigned int n_matches = 0;
 
     // Match with a complexity of O(n^2). This could be improved...
-    std::bitset<3> orientation;
     using PairIterator =
       typename std::set<std::pair<CellIterator, unsigned int>>::const_iterator;
     for (PairIterator it1 = pairs1.begin(); it1 != pairs1.end(); ++it1)
@@ -2153,18 +2153,21 @@ namespace GridTools
             const CellIterator cell2     = it2->first;
             const unsigned int face_idx1 = it1->second;
             const unsigned int face_idx2 = it2->second;
-            if (GridTools::orthogonal_equality(orientation,
-                                               cell1->face(face_idx1),
-                                               cell2->face(face_idx2),
-                                               direction,
-                                               offset,
-                                               matrix))
+            if (const std::optional<std::bitset<3>> orientation =
+                  GridTools::orthogonal_equality(cell1->face(face_idx1),
+                                                 cell2->face(face_idx2),
+                                                 direction,
+                                                 offset,
+                                                 matrix))
               {
                 // We have a match, so insert the matching pairs and
                 // remove the matched cell in pairs2 to speed up the
                 // matching:
                 const PeriodicFacePair<CellIterator> matched_face = {
-                  {cell1, cell2}, {face_idx1, face_idx2}, orientation, matrix};
+                  {cell1, cell2},
+                  {face_idx1, face_idx2},
+                  orientation.value(),
+                  matrix};
                 matched_pairs.push_back(matched_face);
                 pairs2.erase(it2);
                 ++n_matches;
@@ -2504,9 +2507,8 @@ namespace GridTools
 
 
   template <typename FaceIterator>
-  inline bool
+  std::optional<std::bitset<3>>
   orthogonal_equality(
-    std::bitset<3>                                               &orientation,
     const FaceIterator                                           &face1,
     const FaceIterator                                           &face2,
     const unsigned int                                            direction,
@@ -2548,15 +2550,10 @@ namespace GridTools
       }
 
     // And finally, a lookup to determine the ordering bitmask:
-    if (face2_vertices.empty())
-      orientation = OrientationLookupTable<dim>::lookup(matching);
-
-    return face2_vertices.empty();
-
+    return face2_vertices.empty() ?
+             std::make_optional(OrientationLookupTable<dim>::lookup(matching)) :
+             std::nullopt;
   }
-
-
-
 } // namespace GridTools
 
 
index d08016bad84457c9ff75b3d18a7c7535f3f6632e..c5fb404c69d8c7477438fc19bef60896a04685b8 100644 (file)
@@ -250,18 +250,16 @@ for (X : SEQUENTIAL_TRIANGULATION_AND_DOFHANDLER;
     namespace GridTools
     \{
 
-      template bool
+      template std::optional<std::bitset<3>>
       orthogonal_equality<X::active_face_iterator>(
-        std::bitset<3> &,
         const X::active_face_iterator &,
         const X::active_face_iterator &,
         const unsigned int,
         const Tensor<1, deal_II_space_dimension> &,
         const FullMatrix<double> &);
 
-      template bool
+      template std::optional<std::bitset<3>>
       orthogonal_equality<X::face_iterator>(
-        std::bitset<3> &,
         const X::face_iterator &,
         const X::face_iterator &,
         const unsigned int,
index 0e42ae238cb25b0540bfdefe6161fe7cd0562bd4..adc85f55021972f076f8920aa7b42cd6be2e2f38 100644 (file)
@@ -317,24 +317,14 @@ print_matching(DoFHandler<dim, spacedim> &dof_handler,
   deallog << std::endl;
 
 
-  std::bitset<3> orientation;
-  if (not GridTools::orthogonal_equality(orientation,
-                                         face_1,
-                                         face_2,
-                                         dim == 2 ? 1 : 2,
-                                         dealii::Tensor<1, spacedim>()))
-    std::cerr << " not match! oh noze!! " << std::endl;
-  deallog << "Orientation: " << orientation[0] << orientation[1]
-          << orientation[2] << std::endl;
-
-
-  DoFTools::make_periodicity_constraints(face_1,
-                                         face_2,
-                                         constraint_matrix,
-                                         velocity_mask,
-                                         orientation[0],
-                                         orientation[1],
-                                         orientation[2]);
+  const auto orientation = GridTools::orthogonal_equality(
+    face_1, face_2, dim == 2 ? 1 : 2, dealii::Tensor<1, spacedim>());
+  AssertThrow(orientation, ExcMessage(" not match! oh noze!! "));
+  const auto o = *orientation;
+  deallog << "Orientation: " << o[0] << o[1] << o[2] << std::endl;
+
+  DoFTools::make_periodicity_constraints(
+    face_1, face_2, constraint_matrix, velocity_mask, o[0], o[1], o[2]);
   deallog << "Matching:" << std::endl;
   constraint_matrix.print(deallog.get_file_stream());
   constraint_matrix.close();
@@ -343,11 +333,9 @@ print_matching(DoFHandler<dim, spacedim> &dof_handler,
                                          face_1,
                                          constraint_matrix_reverse,
                                          velocity_mask,
-                                         orientation[0],
-                                         orientation[0] ?
-                                           orientation[1] ^ orientation[2] :
-                                           orientation[1],
-                                         orientation[2]);
+                                         o[0],
+                                         o[0] ? o[1] ^ o[2] : o[1],
+                                         o[2]);
   deallog << "Reverse Matching:" << std::endl;
   constraint_matrix_reverse.print(deallog.get_file_stream());
   constraint_matrix_reverse.close();
index 6b2684124778667fda0d09e8a22cd4e5ce3fb0cb..116796598b4320c3bef127b7660991fbf7db533b 100644 (file)
@@ -323,24 +323,14 @@ print_matching(DoFHandler<dim, spacedim> &dof_handler,
   deallog << std::endl;
 
 
-  std::bitset<3> orientation;
-  if (not GridTools::orthogonal_equality(orientation,
-                                         face_1,
-                                         face_2,
-                                         dim == 2 ? 1 : 2,
-                                         dealii::Tensor<1, spacedim>()))
-    std::cerr << " not match! oh noze!! " << std::endl;
-  deallog << "Orientation: " << orientation[0] << orientation[1]
-          << orientation[2] << std::endl;
-
-
-  DoFTools::make_periodicity_constraints(face_1,
-                                         face_2,
-                                         constraint_matrix,
-                                         velocity_mask,
-                                         orientation[0],
-                                         orientation[1],
-                                         orientation[2]);
+  const auto orientation = GridTools::orthogonal_equality(
+    face_1, face_2, dim == 2 ? 1 : 2, dealii::Tensor<1, spacedim>());
+  AssertThrow(orientation, ExcMessage(" not match! oh noze!! "));
+  const auto o = *orientation;
+  deallog << "Orientation: " << o[0] << o[1] << o[2] << std::endl;
+
+  DoFTools::make_periodicity_constraints(
+    face_1, face_2, constraint_matrix, velocity_mask, o[0], o[1], o[2]);
   deallog << "Matching:" << std::endl;
   constraint_matrix.print(deallog.get_file_stream());
   constraint_matrix.close();
@@ -349,11 +339,9 @@ print_matching(DoFHandler<dim, spacedim> &dof_handler,
                                          face_1,
                                          constraint_matrix_reverse,
                                          velocity_mask,
-                                         orientation[0],
-                                         orientation[0] ?
-                                           orientation[1] ^ orientation[2] :
-                                           orientation[1],
-                                         orientation[2]);
+                                         o[0],
+                                         o[0] ? o[1] ^ o[2] : o[1],
+                                         o[2]);
   deallog << "Reverse Matching:" << std::endl;
   constraint_matrix_reverse.print(deallog.get_file_stream());
   constraint_matrix_reverse.close();

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.