From: Johannes Heinz Date: Mon, 4 Jul 2022 07:31:06 +0000 (+0200) Subject: copy boundary ids to remaining boundary faces X-Git-Tag: v9.5.0-rc1~1097^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=713608d06939e484feefc41f451d6e238791f891;p=dealii.git copy boundary ids to remaining boundary faces Co-authored-by: Peter Munch --- diff --git a/doc/news/changes/minor/20220707Heinz b/doc/news/changes/minor/20220707Heinz new file mode 100644 index 0000000000..26dfdb4db9 --- /dev/null +++ b/doc/news/changes/minor/20220707Heinz @@ -0,0 +1,3 @@ +Improved: The function GridGenerator::merge_triangulations() can now copy boundary ids to remaining boundary faces. +
+(Johannes Heinz, 2022/07/07) diff --git a/include/deal.II/grid/grid_generator.h b/include/deal.II/grid/grid_generator.h index d6de5006ab..73865b5f29 100644 --- a/include/deal.II/grid/grid_generator.h +++ b/include/deal.II/grid/grid_generator.h @@ -1884,10 +1884,8 @@ namespace GridGenerator * * @note The function copies the material ids of the cells of the two input * triangulations into the output triangulation. If @p copy_manifold_ids is - * set to @p true, manifold ids will be copied. Boundary indicators are never - * copied. In other words, if the two coarse meshes have anything but the - * default boundary indicators, then you will have to set boundary indicators - * again by hand in the output triangulation. + * set to @p true, manifold ids will be copied. If @p copy_boundary_ids is + * set to @p true, boundary_ids are copied to all remaining faces at the boundary. * * @note This function does not attach any manifolds to @p result, nor does * it set any manifold ids. In particular, manifolds attached to the two @@ -1903,7 +1901,8 @@ namespace GridGenerator const Triangulation &triangulation_2, Triangulation & result, const double duplicated_vertex_tolerance = 1.0e-12, - const bool copy_manifold_ids = false); + const bool copy_manifold_ids = false, + const bool copy_boundary_ids = false); /** * Same as above but allows to merge more than two triangulations at once. @@ -1916,6 +1915,7 @@ namespace GridGenerator * GridGenerator::merge_triangulations({&tria_1, &tria_2, &tria_3}, * merged_triangulation, * 1.0e-10, + * false, * false); * @endcode */ @@ -1925,7 +1925,8 @@ namespace GridGenerator const std::vector *> &triangulations, Triangulation & result, const double duplicated_vertex_tolerance = 1.0e-12, - const bool copy_manifold_ids = false); + const bool copy_manifold_ids = false, + const bool copy_boundary_ids = false); /** * \brief Replicate a given triangulation in multiple coordinate axes diff --git a/source/grid/grid_generator.cc b/source/grid/grid_generator.cc index 41fa3a00f4..637a714d0e 100644 --- a/source/grid/grid_generator.cc +++ b/source/grid/grid_generator.cc @@ -6469,7 +6469,8 @@ namespace GridGenerator const std::vector *> &triangulations, Triangulation & result, const double duplicated_vertex_tolerance, - const bool copy_manifold_ids) + const bool copy_manifold_ids, + const bool copy_boundary_ids) { std::vector> vertices; std::vector> cells; @@ -6545,6 +6546,22 @@ namespace GridGenerator result.clear(); result.create_triangulation(vertices, cells, subcell_data); + if (copy_boundary_ids) + { + auto result_cell = result.begin(); + for (const auto &tria : triangulations) + { + for (const auto &cell : tria->cell_iterators()) + { + for (auto const &f : cell->face_indices()) + if (result_cell->face(f)->at_boundary()) + result_cell->face(f)->set_boundary_id( + cell->face(f)->boundary_id()); + result_cell++; + } + } + } + Assert(duplicated_vertex_tolerance > 0.0 || n_accumulated_vertices == result.n_vertices(), ExcInternalError()); @@ -6558,7 +6575,8 @@ namespace GridGenerator const Triangulation &triangulation_2, Triangulation & result, const double duplicated_vertex_tolerance, - const bool copy_manifold_ids) + const bool copy_manifold_ids, + const bool copy_boundary_ids) { // if either Triangulation is empty then merging is just a copy. if (triangulation_1.n_cells() == 0) @@ -6574,7 +6592,8 @@ namespace GridGenerator merge_triangulations({&triangulation_1, &triangulation_2}, result, duplicated_vertex_tolerance, - copy_manifold_ids); + copy_manifold_ids, + copy_boundary_ids); } diff --git a/source/grid/grid_generator.inst.in b/source/grid/grid_generator.inst.in index 05d1fff1c2..c0b6cce0bd 100644 --- a/source/grid/grid_generator.inst.in +++ b/source/grid/grid_generator.inst.in @@ -84,6 +84,7 @@ for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS) const Triangulation *> &, Triangulation &, const double, + const bool, const bool); template void @@ -94,7 +95,8 @@ for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS) &triangulation_2, Triangulation &result, const double duplicated_vertex_tolerance, - const bool copy_manifold_ids); + const bool copy_manifold_ids, + const bool copy_boundary_ids); template void replicate_triangulation( diff --git a/tests/grid/grid_tools_07.cc b/tests/grid/grid_tools_07.cc index a40c3c3667..141a680ed5 100644 --- a/tests/grid/grid_tools_07.cc +++ b/tests/grid/grid_tools_07.cc @@ -65,21 +65,7 @@ generate_grid(Triangulation<2> &triangulation) // tolerance 0. to ensure vertices are not merged GridGenerator::merge_triangulations( - triangulation1, triangulation2, triangulation, 0., false); - - // make sure boundary ids are kept - std::vector boundary_ids; - for (const auto &face : triangulation1.active_face_iterators()) - if (face->at_boundary()) - boundary_ids.emplace_back(face->boundary_id()); - for (const auto &face : triangulation2.active_face_iterators()) - if (face->at_boundary()) - boundary_ids.emplace_back(face->boundary_id()); - - unsigned int i = 0; - for (const auto &face : triangulation.active_face_iterators()) - if (face->at_boundary()) - face->set_boundary_id(boundary_ids[i++]); + triangulation1, triangulation2, triangulation, 0., false, true); } std::vector diff --git a/tests/grid/merge_triangulations_10.cc b/tests/grid/merge_triangulations_10.cc new file mode 100644 index 0000000000..d7011278ae --- /dev/null +++ b/tests/grid/merge_triangulations_10.cc @@ -0,0 +1,73 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2018 - 2020 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 for flag copy_boundary_ids in merge_triangulations() + +#include +#include + +#include "../tests.h" + +template +void +print_boundary_ids(const dealii::Triangulation &tria) +{ + for (const auto &face : tria.active_face_iterators()) + if (face->at_boundary()) + deallog << "Face at " << face->center() + << " has id: " << face->boundary_id() << "\n"; + deallog << "-----------------" << std::endl; +} + + +int +main() +{ + using namespace dealii; + + initlog(); + + const unsigned int dim = 2; + + Triangulation tria1; + Triangulation tria2; + Triangulation tria3; + Triangulation tria; + GridGenerator::hyper_rectangle(tria1, {-1., -1.}, {0., 1.}, true); + GridGenerator::hyper_rectangle(tria2, {0., -1.}, {1., 1.}, true); + for (auto &face : tria2.active_face_iterators()) + if (face->at_boundary()) + face->set_boundary_id(face->boundary_id() + dim * 2); + + deallog << "Initial IDs of Triangulation 1" << std::endl; + print_boundary_ids(tria1); + deallog << "Initial IDs of Triangulation 2" << std::endl; + print_boundary_ids(tria2); + + deallog << "Merge two rectangles, one face is eliminated" << std::endl; + GridGenerator::merge_triangulations(tria1, tria2, tria, 1e-6, false, true); + print_boundary_ids(tria); + + deallog + << "Merge two rectangles, due to hanging node there is no internal face" + << std::endl; + tria.clear(); + tria2.refine_global(1); + GridGenerator::flatten_triangulation(tria2, tria3); + GridGenerator::merge_triangulations(tria1, tria3, tria, 1e-6, false, true); + print_boundary_ids(tria); + + deallog << std::endl; +} diff --git a/tests/grid/merge_triangulations_10.output b/tests/grid/merge_triangulations_10.output new file mode 100644 index 0000000000..7c776527c0 --- /dev/null +++ b/tests/grid/merge_triangulations_10.output @@ -0,0 +1,36 @@ + +DEAL::Initial IDs of Triangulation 1 +DEAL::Face at -0.500000 -1.00000 has id: 2 +Face at -1.00000 0.00000 has id: 0 +Face at 0.00000 0.00000 has id: 1 +Face at -0.500000 1.00000 has id: 3 +----------------- +DEAL::Initial IDs of Triangulation 2 +DEAL::Face at 0.500000 -1.00000 has id: 6 +Face at 0.00000 0.00000 has id: 4 +Face at 1.00000 0.00000 has id: 5 +Face at 0.500000 1.00000 has id: 7 +----------------- +DEAL::Merge two rectangles, one face is eliminated +DEAL::Face at -0.500000 -1.00000 has id: 2 +Face at -1.00000 0.00000 has id: 0 +Face at 0.500000 -1.00000 has id: 6 +Face at -0.500000 1.00000 has id: 3 +Face at 0.500000 1.00000 has id: 7 +Face at 1.00000 0.00000 has id: 5 +----------------- +DEAL::Merge two rectangles, due to hanging node there is no internal face +DEAL::Face at -0.500000 -1.00000 has id: 2 +Face at -1.00000 0.00000 has id: 0 +Face at 0.00000 0.00000 has id: 1 +Face at 0.250000 -1.00000 has id: 6 +Face at 0.00000 -0.500000 has id: 4 +Face at -0.500000 1.00000 has id: 3 +Face at 0.250000 1.00000 has id: 7 +Face at 1.00000 -0.500000 has id: 5 +Face at 0.750000 -1.00000 has id: 6 +Face at 0.00000 0.500000 has id: 4 +Face at 1.00000 0.500000 has id: 5 +Face at 0.750000 1.00000 has id: 7 +----------------- +DEAL::