]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Copy periodic faces when using copy_triangulation
authorLaura Prieto Saavedra <lauraprietosaavedra@gmail.com>
Thu, 19 Oct 2023 21:34:00 +0000 (17:34 -0400)
committerLaura Prieto Saavedra <lauraprietosaavedra@gmail.com>
Thu, 19 Oct 2023 21:34:00 +0000 (17:34 -0400)
doc/news/changes/minor/20231019PrietoSaavedraMunch [new file with mode: 0644]
source/grid/tria.cc
tests/grid/copy_triangulation_pbc.cc [new file with mode: 0644]
tests/grid/copy_triangulation_pbc.output [new file with mode: 0644]

diff --git a/doc/news/changes/minor/20231019PrietoSaavedraMunch b/doc/news/changes/minor/20231019PrietoSaavedraMunch
new file mode 100644 (file)
index 0000000..1d83009
--- /dev/null
@@ -0,0 +1,4 @@
+Fixed: The periodic faces information was not being copied
+when using the `Triangulation::copy_triangulation()` function. This is now fixed.
+<br>
+(Laura Prieto Saavedra, Peter Munch, 2023/10/19)
\ No newline at end of file
index 2ab3f895a4a8a0452d72d00e529cc92f02d7fed6..38c8beed8ca86b8961b07bfee4069cceabe485ec 100644 (file)
@@ -12520,6 +12520,34 @@ void Triangulation<dim, spacedim>::copy_triangulation(
   if (other_tria.policy)
     this->policy = other_tria.policy->clone();
 
+  // periodic faces
+  this->periodic_face_pairs_level_0.reserve(
+    other_tria.periodic_face_pairs_level_0.size());
+
+  for (auto entry : other_tria.periodic_face_pairs_level_0)
+    {
+      entry.cell[0] =
+        cell_iterator(this, entry.cell[0]->level(), entry.cell[0]->index());
+      entry.cell[1] =
+        cell_iterator(this, entry.cell[1]->level(), entry.cell[1]->index());
+      periodic_face_pairs_level_0.emplace_back(entry);
+    }
+
+  for (auto [first_cell_, second_cell_and_orientation] :
+       other_tria.periodic_face_map)
+    {
+      auto first_cell  = first_cell_; // make copy since key is const
+      first_cell.first = cell_iterator(this,
+                                       first_cell.first->level(),
+                                       first_cell.first->index());
+      second_cell_and_orientation.first.first =
+        cell_iterator(this,
+                      second_cell_and_orientation.first.first->level(),
+                      second_cell_and_orientation.first.first->index());
+
+      this->periodic_face_map[first_cell] = second_cell_and_orientation;
+    }
+
   // inform those who are listening on other_tria of the copy operation
   other_tria.signals.copy(*this);
   // also inform all listeners of the current triangulation that the
diff --git a/tests/grid/copy_triangulation_pbc.cc b/tests/grid/copy_triangulation_pbc.cc
new file mode 100644 (file)
index 0000000..f31a791
--- /dev/null
@@ -0,0 +1,75 @@
+/* ---------------------------------------------------------------------
+ *
+ * Copyright (C) 2023 by the deal.II authors
+ *
+ * This file is part of the deal.II library.
+ *
+ * The deal.II library is free software; you can use it, redistribute
+ * it, and/or modify it under the terms of the GNU Lesser General
+ * Public License as published by the Free Software Foundation; either
+ * version 2.1 of the License, or (at your option) any later version.
+ * The full text of the license can be found in the file LICENSE.md at
+ * the top level directory of deal.II.
+ *
+ * ---------------------------------------------------------------------
+ *
+ * Test that when copy_triangulation is used also the information about
+ * periodic faces is copied.
+ */
+#include <deal.II/dofs/dof_handler.h>
+#include <deal.II/dofs/dof_tools.h>
+
+#include <deal.II/fe/fe_q.h>
+#include <deal.II/fe/fe_system.h>
+
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/grid_tools.h>
+#include <deal.II/grid/tria.h>
+
+#include "../tests.h"
+
+int
+main()
+{
+  initlog();
+
+  const unsigned int dim           = 2;
+  const unsigned int n_components  = 1;
+  const unsigned int n_refinements = 1;
+
+  Triangulation<dim> tria;
+  GridGenerator::hyper_cube(tria, 0, 1, true);
+
+  std::vector<
+    GridTools::PeriodicFacePair<typename Triangulation<dim>::cell_iterator>>
+    periodic_faces;
+
+  GridTools::collect_periodic_faces(tria, 0, 1, 0, periodic_faces);
+
+  tria.add_periodicity(periodic_faces);
+  tria.refine_global(n_refinements);
+
+  DoFHandler<dim> dof_handler(tria);
+  dof_handler.distribute_dofs(FESystem<dim>(FE_Q<dim>(1), n_components));
+
+  deallog << "Tria:" << std::endl;
+  for (const auto &elem : tria.get_periodic_face_map())
+    {
+      deallog << elem.first.first << " " << elem.first.second << " "
+              << elem.second.first.first << " " << elem.second.first.second
+              << std::endl;
+    }
+
+  Triangulation<dim> other_tria;
+  other_tria.copy_triangulation(tria);
+
+  deallog << "Copy of tria:" << std::endl;
+  for (const auto &elem : other_tria.get_periodic_face_map())
+    {
+      deallog << elem.first.first << " " << elem.first.second << " "
+              << elem.second.first.first << " " << elem.second.first.second
+              << std::endl;
+    }
+
+  deallog << "OK" << std::endl;
+}
diff --git a/tests/grid/copy_triangulation_pbc.output b/tests/grid/copy_triangulation_pbc.output
new file mode 100644 (file)
index 0000000..aa0778c
--- /dev/null
@@ -0,0 +1,16 @@
+
+DEAL::Tria:
+DEAL::0.0 0 0.0 1 
+DEAL::0.0 1 0.0 0 
+DEAL::1.0 0 1.1 1 
+DEAL::1.1 1 1.0 0 
+DEAL::1.2 0 1.3 1 
+DEAL::1.3 1 1.2 0 
+DEAL::Copy of tria:
+DEAL::0.0 0 0.0 1 
+DEAL::0.0 1 0.0 0 
+DEAL::1.0 0 1.1 1 
+DEAL::1.1 1 1.0 0 
+DEAL::1.2 0 1.3 1 
+DEAL::1.3 1 1.2 0
+DEAL::OK

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.