From: Laura Prieto Saavedra Date: Thu, 19 Oct 2023 21:34:00 +0000 (-0400) Subject: Copy periodic faces when using copy_triangulation X-Git-Tag: relicensing~376^2~1 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=f0f920b2c66ab33e61ce3857477b181e82f661cb;p=dealii.git Copy periodic faces when using copy_triangulation --- diff --git a/doc/news/changes/minor/20231019PrietoSaavedraMunch b/doc/news/changes/minor/20231019PrietoSaavedraMunch new file mode 100644 index 0000000000..1d830099aa --- /dev/null +++ b/doc/news/changes/minor/20231019PrietoSaavedraMunch @@ -0,0 +1,4 @@ +Fixed: The periodic faces information was not being copied +when using the `Triangulation::copy_triangulation()` function. This is now fixed. +
+(Laura Prieto Saavedra, Peter Munch, 2023/10/19) \ No newline at end of file diff --git a/source/grid/tria.cc b/source/grid/tria.cc index 2ab3f895a4..38c8beed8c 100644 --- a/source/grid/tria.cc +++ b/source/grid/tria.cc @@ -12520,6 +12520,34 @@ void Triangulation::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 index 0000000000..f31a791479 --- /dev/null +++ b/tests/grid/copy_triangulation_pbc.cc @@ -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 +#include + +#include +#include + +#include +#include +#include + +#include "../tests.h" + +int +main() +{ + initlog(); + + const unsigned int dim = 2; + const unsigned int n_components = 1; + const unsigned int n_refinements = 1; + + Triangulation tria; + GridGenerator::hyper_cube(tria, 0, 1, true); + + std::vector< + GridTools::PeriodicFacePair::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 dof_handler(tria); + dof_handler.distribute_dofs(FESystem(FE_Q(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 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 index 0000000000..aa0778c999 --- /dev/null +++ b/tests/grid/copy_triangulation_pbc.output @@ -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