]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Remove deprecated Triangulation::create_triangulation_compatibility
authorDaniel Arndt <arndtd@ornl.gov>
Sun, 2 Jul 2023 13:14:49 +0000 (09:14 -0400)
committerDaniel Arndt <arndtd@ornl.gov>
Mon, 3 Jul 2023 21:51:15 +0000 (17:51 -0400)
13 files changed:
examples/step-1/doc/intro.dox
include/deal.II/grid/persistent_tria.h
include/deal.II/grid/tria.h
source/grid/persistent_tria.cc
source/grid/tria.cc
tests/bits/neighboring_cells_at_two_faces.cc
tests/bits/q_points.cc
tests/data_out/data_out_curved_cells.cc
tests/grid/grid_invert_02.cc
tests/grid/line_coarsening_3d.cc
tests/grid/mesh_3d.h
tests/grid/subcelldata.cc
tests/lac/constraints.cc

index ead99932b8021002aecc69e223240c52c1835d58..c0a2b20bae7b3c68cee140bbdc1c9442b390218f 100644 (file)
@@ -17,8 +17,8 @@ different levels:
   every single (member) function in deal.II. You get there if, for
   example, you click on the "Main page" or "Classes" tab at the top of
   this page. This is the place where you would look up what the second
-  argument of Triangulation::create_triangulation_compatibility means,
-  to give just one slightly obscure example. You need this level of
+  argument of Triangulation::create_triangulation means,
+  to give just one example. You need this level of
   documentation for when you know what you want to do, but forgot how
   exactly the function was named, what its arguments are, or what it
   returns. Note that you also get into the manual whenever you read
index 8c99d4951447cba24ae9be2c8f18cfe333cb9bfe..7ee795b52e236e1a5fbefe444da739fc4717e0c0 100644 (file)
@@ -211,19 +211,6 @@ public:
     const TriangulationDescription::Description<dim, spacedim>
       &construction_data) override;
 
-  /**
-   * An overload of the respective function of the base class.
-   *
-   * Throw an error, since this function is not useful in the context of this
-   * class.
-   */
-  DEAL_II_DEPRECATED
-  virtual void
-  create_triangulation_compatibility(
-    const std::vector<Point<spacedim>> &vertices,
-    const std::vector<CellData<dim>> &  cells,
-    const SubCellData &                 subcelldata) override;
-
   /**
    * Write all refine and coarsen flags to the ostream @p out.
    */
index f3b01b39767cbb8a18cfef2a17c3e8e04bd54230..726e3bcad707d1ca737eda25ded8115ddbecc5d8 100644 (file)
@@ -1893,21 +1893,6 @@ public:
     const TriangulationDescription::Description<dim, spacedim>
       &construction_data);
 
-  /**
-   * For backward compatibility, only. This function takes the cell data in
-   * the ordering as requested by deal.II versions up to 5.2, converts it to
-   * the new (lexicographic) ordering and calls create_triangulation().
-   *
-   * @note This function internally calls create_triangulation and therefore
-   * can throw the same exception as the other function.
-   */
-  DEAL_II_DEPRECATED
-  virtual void
-  create_triangulation_compatibility(
-    const std::vector<Point<spacedim>> &vertices,
-    const std::vector<CellData<dim>> &  cells,
-    const SubCellData &                 subcelldata);
-
   /**
    * Revert or flip the direction_flags of a dim<spacedim triangulation, see
    * @ref GlossDirectionFlag.
index 21189ab35f3c2fedea2ed18d4e79fe6b10fe6ffa..a5913ee39e44c97853dcab5d1ed0432d46ca75a7 100644 (file)
@@ -157,18 +157,6 @@ PersistentTriangulation<dim, spacedim>::create_triangulation(
 
 
 
-template <int dim, int spacedim>
-void
-PersistentTriangulation<dim, spacedim>::create_triangulation_compatibility(
-  const std::vector<Point<spacedim>> &,
-  const std::vector<CellData<dim>> &,
-  const SubCellData &)
-{
-  Assert(false, ExcImpossibleInDim(dim));
-}
-
-
-
 template <int dim, int spacedim>
 void
 PersistentTriangulation<dim, spacedim>::write_flags(std::ostream &out) const
index c5682522294637f8a8594e538d05afc9454152fe..73570ab6c3596abed0caae48580c3662208b64f6 100644 (file)
@@ -414,59 +414,6 @@ namespace
 
 
 
-  /**
-   * A set of three functions that
-   * reorder the data given to
-   * create_triangulation_compatibility
-   * from the "classic" to the
-   * "current" format of vertex
-   * numbering of cells and
-   * faces. These functions do the
-   * reordering of their arguments
-   * in-place.
-   */
-  void
-  reorder_compatibility(const std::vector<CellData<1>> &, const SubCellData &)
-  {
-    // nothing to do here: the format
-    // hasn't changed for 1d
-  }
-
-
-  void
-  reorder_compatibility(std::vector<CellData<2>> &cells, const SubCellData &)
-  {
-    for (auto &cell : cells)
-      if (cell.vertices.size() == GeometryInfo<2>::vertices_per_cell)
-        std::swap(cell.vertices[2], cell.vertices[3]);
-  }
-
-
-  void
-  reorder_compatibility(std::vector<CellData<3>> &cells,
-                        SubCellData &             subcelldata)
-  {
-    unsigned int tmp[GeometryInfo<3>::vertices_per_cell];
-    static constexpr std::array<unsigned int,
-                                GeometryInfo<3>::vertices_per_cell>
-      local_vertex_numbering{{0, 1, 5, 4, 2, 3, 7, 6}};
-    for (auto &cell : cells)
-      if (cell.vertices.size() == GeometryInfo<3>::vertices_per_cell)
-        {
-          for (const unsigned int i : GeometryInfo<3>::vertex_indices())
-            tmp[i] = cell.vertices[i];
-          for (const unsigned int i : GeometryInfo<3>::vertex_indices())
-            cell.vertices[local_vertex_numbering[i]] = tmp[i];
-        }
-
-    // now points in boundary quads
-    for (auto &boundary_quad : subcelldata.boundary_quads)
-      if (boundary_quad.vertices.size() == GeometryInfo<2>::vertices_per_cell)
-        std::swap(boundary_quad.vertices[2], boundary_quad.vertices[3]);
-  }
-
-
-
   /**
    * Return the index of the vertex
    * in the middle of this object,
@@ -11548,25 +11495,6 @@ void Triangulation<dim, spacedim>::copy_triangulation(
 
 
 
-template <int dim, int spacedim>
-DEAL_II_CXX20_REQUIRES((concepts::is_valid_dim_spacedim<dim, spacedim>))
-void Triangulation<dim, spacedim>::create_triangulation_compatibility(
-  const std::vector<Point<spacedim>> &v,
-  const std::vector<CellData<dim>> &  cells,
-  const SubCellData &                 subcelldata)
-{
-  std::vector<CellData<dim>> reordered_cells(cells);             // NOLINT
-  SubCellData                reordered_subcelldata(subcelldata); // NOLINT
-
-  // in-place reordering of data
-  reorder_compatibility(reordered_cells, reordered_subcelldata);
-
-  // now create triangulation from
-  // reordered data
-  create_triangulation(v, reordered_cells, reordered_subcelldata);
-}
-
-
 template <int dim, int spacedim>
 DEAL_II_CXX20_REQUIRES((concepts::is_valid_dim_spacedim<dim, spacedim>))
 void Triangulation<dim, spacedim>::reset_policy()
index 4ad6151b44f7f5e91196398bc72d2288dbd1857e..b675f4aa30d8903a37bb80a4d604eb108253f1d0 100644 (file)
@@ -60,19 +60,19 @@ create_grid(Triangulation<2> &tria)
   std::vector<CellData<2>> cells(2);
   cells[0].vertices[0] = 0;
   cells[0].vertices[1] = 1;
-  cells[0].vertices[2] = 2;
-  cells[0].vertices[3] = 3;
+  cells[0].vertices[2] = 3;
+  cells[0].vertices[3] = 2;
   cells[1].vertices[0] = 0;
   cells[1].vertices[1] = 4;
-  cells[1].vertices[2] = 2;
-  cells[1].vertices[3] = 1;
+  cells[1].vertices[2] = 1;
+  cells[1].vertices[3] = 2;
 
   // generate a triangulation
   // out of this
   GridTools::consistently_order_cells(cells);
   try
     {
-      tria.create_triangulation_compatibility(vertices, cells, SubCellData());
+      tria.create_triangulation(vertices, cells, SubCellData());
     }
   catch (Triangulation<2>::DistortedCellList &dcv)
     {
index d5b04dd97a9dd144e1149bc281f21d66e8c4c16f..377856ea096286e8b19fcaa4c2e5008ba6b4069e 100644 (file)
 
 #include "../tests.h"
 
-
+namespace
+{
+  /**
+   * This file uses a different ordering for the vertices in a hex
+   * cell than we usually do in deal.II. The different convention used
+   * here originates in what we believed the ordering to be in UCD
+   * format, until it was discovered in 2022 that UCD will interpret
+   * this ordering to correspond to inverted cells -- as a
+   * consequence, the UCD ordering was fixed, but the current file is
+   * stuck on the old ordering.
+   */
+  constexpr std::array<unsigned int, 8> local_vertex_numbering{
+    {0, 1, 5, 4, 2, 3, 7, 6}};
+
+  /**
+   * Following is a set of functions that reorder the data from the
+   * "current" to the "classic" format of vertex numbering of cells
+   * and faces. These functions do the reordering of their arguments
+   * in-place.
+   */
+  void
+  reorder_old_to_new_style(std::vector<CellData<2>> &cells)
+  {
+    for (auto &cell : cells)
+      std::swap(cell.vertices[2], cell.vertices[3]);
+  }
+
+
+  void
+  reorder_old_to_new_style(std::vector<CellData<3>> &cells)
+  {
+    // undo the ordering above
+    unsigned int tmp[GeometryInfo<3>::vertices_per_cell];
+    for (auto &cell : cells)
+      {
+        for (const unsigned int i : GeometryInfo<3>::vertex_indices())
+          tmp[i] = cell.vertices[i];
+        for (const unsigned int i : GeometryInfo<3>::vertex_indices())
+          cell.vertices[local_vertex_numbering[i]] = tmp[i];
+      }
+  }
+} // namespace
 
 void
 create_two_cubes(Triangulation<3> &coarse_grid)
@@ -69,10 +110,9 @@ create_two_cubes(Triangulation<3> &coarse_grid)
 
   // finally generate a triangulation
   // out of this
+  reorder_old_to_new_style(cells);
   GridTools::consistently_order_cells(cells);
-  coarse_grid.create_triangulation_compatibility(vertices,
-                                                 cells,
-                                                 SubCellData());
+  coarse_grid.create_triangulation(vertices, cells, SubCellData());
 }
 
 
index cea603b2999ec2bf263823d86b986f6f181404d6..de957fb0920b94ca3e6eb576d290052e2b2a77a6 100644 (file)
@@ -124,14 +124,12 @@ curved_grid(std::ostream &out)
 
         cells[index].vertices[0] = p_iplus;
         cells[index].vertices[1] = p;
-        cells[index].vertices[2] = p_jplus;
-        cells[index].vertices[3] = p_iplus_jplus;
+        cells[index].vertices[2] = p_iplus_jplus;
+        cells[index].vertices[3] = p_jplus;
       }
   // create triangulation
   Triangulation<2> triangulation;
-  triangulation.create_triangulation_compatibility(vertices,
-                                                   cells,
-                                                   SubCellData());
+  triangulation.create_triangulation(vertices, cells, SubCellData());
   // now provide everything that is
   // needed for solving a Laplace
   // equation.
index 5cf5a6bc5bb2774dc750a219d76f94ccb116f678..0a3ce5743feffa6133c3d302ccce7ba7657ffcf4 100644 (file)
@@ -58,8 +58,11 @@ test()
 
   deallog << "We inverted " << n_cells_inverted << " cells." << std::endl;
 
+  std::swap(cells[0].vertices[2], cells[0].vertices[3]);
+  std::swap(cells[1].vertices[2], cells[1].vertices[3]);
+
   Triangulation<dim> tria;
-  tria.create_triangulation_compatibility(vertices, cells, subcelldata);
+  tria.create_triangulation(vertices, cells, subcelldata);
 
   std::ostream &logfile = deallog.get_file_stream();
   logfile << "---------------------------------------------" << std::endl
index ed7994d831df298e7529bc07b3b1dbd714e6de70..20a12d6bc6a9bbf2c18a3a1fc3922cd4a278f28a 100644 (file)
@@ -55,16 +55,16 @@ create_star_structured_cylinder(Triangulation<3> & coarse_grid,
     {
       cells[c].vertices[0] = 0;
       cells[c].vertices[1] = 1 + 2 * c;
-      cells[c].vertices[2] = 2 + 2 * c;
-      cells[c].vertices[3] = (3 + 2 * c) % (2 * n_cells);
-      cells[c].vertices[4] = 0 + 2 * n_cells + 1;
-      cells[c].vertices[5] = 1 + 2 * c + 2 * n_cells + 1;
-      cells[c].vertices[6] = 2 + 2 * c + 2 * n_cells + 1;
-      cells[c].vertices[7] = (3 + 2 * c) % (2 * n_cells) + 2 * n_cells + 1;
+      cells[c].vertices[2] = 0 + 2 * n_cells + 1;
+      cells[c].vertices[3] = 1 + 2 * c + 2 * n_cells + 1;
+      cells[c].vertices[4] = (3 + 2 * c) % (2 * n_cells);
+      cells[c].vertices[5] = 2 + 2 * c;
+      cells[c].vertices[6] = (3 + 2 * c) % (2 * n_cells) + 2 * n_cells + 1;
+      cells[c].vertices[7] = 2 + 2 * c + 2 * n_cells + 1;
     }
   // finally generate a triangulation
   // out of this
-  coarse_grid.create_triangulation_compatibility(points, cells, SubCellData());
+  coarse_grid.create_triangulation(points, cells, SubCellData());
 }
 
 
index a7d5d3b317ecd3038f939e9b578644630c2420b1..f4564b3c8ad0aad351c3405c628c1e2defe6ca50 100644 (file)
@@ -139,9 +139,7 @@ create_two_cubes_rotation(Triangulation<3> & coarse_grid,
     }
   // finally generate a triangulation
   // out of this
-  coarse_grid.create_triangulation_compatibility(vertices,
-                                                 cells,
-                                                 SubCellData());
+  coarse_grid.create_triangulation(vertices, cells, SubCellData());
 }
 
 
index 23dce0cc9f833005eba485627ace3caddad299b4..cc26ed861e59213a598a221701b07629f7ff93cc 100644 (file)
 #include "../tests.h"
 
 
+namespace
+{
+  /**
+   * This file uses a different ordering for the vertices in a hex
+   * cell than we usually do in deal.II. The different convention used
+   * here originates in what we believed the ordering to be in UCD
+   * format, until it was discovered in 2022 that UCD will interpret
+   * this ordering to correspond to inverted cells -- as a
+   * consequence, the UCD ordering was fixed, but the current file is
+   * stuck on the old ordering.
+   */
+  constexpr std::array<unsigned int, 8> local_vertex_numbering{
+    {0, 1, 5, 4, 2, 3, 7, 6}};
+
+  /**
+   * And now also in the opposite direction.
+   */
+  void
+  reorder_old_to_new_style(std::vector<CellData<3>> &cells)
+  {
+    // undo the ordering above
+    unsigned int tmp[GeometryInfo<3>::vertices_per_cell];
+    for (auto &cell : cells)
+      {
+        for (const unsigned int i : GeometryInfo<3>::vertex_indices())
+          tmp[i] = cell.vertices[i];
+        for (const unsigned int i : GeometryInfo<3>::vertex_indices())
+          cell.vertices[local_vertex_numbering[i]] = tmp[i];
+      }
+  }
+
+  /**
+   * Following is a set of functions that reorder the data from the
+   * "current" to the "classic" format of vertex numbering of cells
+   * and faces. These functions do the reordering of their arguments
+   * in-place.
+   */
+  void
+  reorder_old_to_new_style(std::vector<CellData<2>> &cells)
+  {
+    for (auto &cell : cells)
+      std::swap(cell.vertices[2], cell.vertices[3]);
+  }
+} // namespace
+
 
 static unsigned subcells[6][4] = {{0, 1, 2, 3},
                                   {4, 5, 6, 7},
@@ -93,8 +138,9 @@ test()
         }
     }
 
+  reorder_old_to_new_style(cells);
   Triangulation<dim> tria;
-  tria.create_triangulation_compatibility(vertices, cells, subcelldata);
+  tria.create_triangulation(vertices, cells, subcelldata);
 
   GridOutFlags::Ucd ucd_flags(true, true);
   GridOut           grid_out;
index f342a66bdf67c2fa740a8a881d1adddfc59b46d8..c9324b5b20966d03ac637500290ace989cdd7214 100644 (file)
 #include "../tests.h"
 
 
+namespace
+{
+  /**
+   * This file uses a different ordering for the vertices in a hex
+   * cell than we usually do in deal.II. The different convention used
+   * here originates in what we believed the ordering to be in UCD
+   * format, until it was discovered in 2022 that UCD will interpret
+   * this ordering to correspond to inverted cells -- as a
+   * consequence, the UCD ordering was fixed, but the current file is
+   * stuck on the old ordering.
+   */
+  constexpr std::array<unsigned int, 8> local_vertex_numbering{
+    {0, 1, 5, 4, 2, 3, 7, 6}};
+
+  /**
+   * Following is a set of functions that reorder the data from the
+   * "current" to the "classic" format of vertex numbering of cells
+   * and faces. These functions do the reordering of their arguments
+   * in-place.
+   */
+  void
+  reorder_old_to_new_style(std::vector<CellData<2>> &cells)
+  {
+    for (auto &cell : cells)
+      std::swap(cell.vertices[2], cell.vertices[3]);
+  }
+
+
+  void
+  reorder_old_to_new_style(std::vector<CellData<3>> &cells)
+  {
+    // undo the ordering above
+    unsigned int tmp[GeometryInfo<3>::vertices_per_cell];
+    for (auto &cell : cells)
+      {
+        for (const unsigned int i : GeometryInfo<3>::vertex_indices())
+          tmp[i] = cell.vertices[i];
+        for (const unsigned int i : GeometryInfo<3>::vertex_indices())
+          cell.vertices[local_vertex_numbering[i]] = tmp[i];
+      }
+  }
+} // namespace
+
 
 void
 make_tria(Triangulation<3> &tria, int step)
@@ -72,10 +115,11 @@ make_tria(Triangulation<3> &tria, int step)
           cells[0].material_id = 0;
           cells[1].material_id = 0;
 
-          tria.create_triangulation_compatibility(
-            std::vector<Point<3>>(&vertices[0], &vertices[12]),
-            cells,
-            SubCellData()); // no boundary information
+          reorder_old_to_new_style(cells);
+          tria.create_triangulation(std::vector<Point<3>>(&vertices[0],
+                                                          &vertices[12]),
+                                    cells,
+                                    SubCellData()); // no boundary information
 
           if (step == 0)
             tria.last_active()->set_refine_flag();
@@ -115,10 +159,11 @@ make_tria(Triangulation<3> &tria, int step)
           cells[0].material_id = 0;
           cells[1].material_id = 0;
 
-          tria.create_triangulation_compatibility(
-            std::vector<Point<3>>(&vertices[0], &vertices[12]),
-            cells,
-            SubCellData()); // no boundary information
+          reorder_old_to_new_style(cells);
+          tria.create_triangulation(std::vector<Point<3>>(&vertices[0],
+                                                          &vertices[12]),
+                                    cells,
+                                    SubCellData()); // no boundary information
 
           if (step == 2)
             tria.last_active()->set_refine_flag();
@@ -158,10 +203,11 @@ make_tria(Triangulation<3> &tria, int step)
           cells[0].material_id = 0;
           cells[1].material_id = 0;
 
-          tria.create_triangulation_compatibility(
-            std::vector<Point<3>>(&vertices[0], &vertices[12]),
-            cells,
-            SubCellData()); // no boundary information
+          reorder_old_to_new_style(cells);
+          tria.create_triangulation(std::vector<Point<3>>(&vertices[0],
+                                                          &vertices[12]),
+                                    cells,
+                                    SubCellData()); // no boundary information
 
           if (step == 4)
             tria.last_active()->set_refine_flag();
@@ -215,10 +261,11 @@ make_tria(Triangulation<3> &tria, int step)
           cells[2].material_id = 0;
           cells[3].material_id = 0;
 
-          tria.create_triangulation_compatibility(
-            std::vector<Point<3>>(&vertices[0], &vertices[18]),
-            cells,
-            SubCellData()); // no boundary information
+          reorder_old_to_new_style(cells);
+          tria.create_triangulation(std::vector<Point<3>>(&vertices[0],
+                                                          &vertices[18]),
+                                    cells,
+                                    SubCellData()); // no boundary information
 
           switch (step)
             {

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.