--- /dev/null
+Improved: The function GridGenerator::merge_triangulations() can now copy boundary ids to remaining boundary faces.
+<br>
+(Johannes Heinz, 2022/07/07)
*
* @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
const Triangulation<dim, spacedim> &triangulation_2,
Triangulation<dim, spacedim> & 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.
* GridGenerator::merge_triangulations({&tria_1, &tria_2, &tria_3},
* merged_triangulation,
* 1.0e-10,
+ * false,
* false);
* @endcode
*/
const std::vector<const Triangulation<dim, spacedim> *> &triangulations,
Triangulation<dim, spacedim> & 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
const std::vector<const Triangulation<dim, spacedim> *> &triangulations,
Triangulation<dim, spacedim> & result,
const double duplicated_vertex_tolerance,
- const bool copy_manifold_ids)
+ const bool copy_manifold_ids,
+ const bool copy_boundary_ids)
{
std::vector<Point<spacedim>> vertices;
std::vector<CellData<dim>> cells;
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());
const Triangulation<dim, spacedim> &triangulation_2,
Triangulation<dim, spacedim> & 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)
merge_triangulations({&triangulation_1, &triangulation_2},
result,
duplicated_vertex_tolerance,
- copy_manifold_ids);
+ copy_manifold_ids,
+ copy_boundary_ids);
}
const Triangulation<deal_II_dimension, deal_II_space_dimension> *> &,
Triangulation<deal_II_dimension, deal_II_space_dimension> &,
const double,
+ const bool,
const bool);
template void
&triangulation_2,
Triangulation<deal_II_dimension, deal_II_space_dimension> &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(
// tolerance 0. to ensure vertices are not merged
GridGenerator::merge_triangulations(
- triangulation1, triangulation2, triangulation, 0., false);
-
- // make sure boundary ids are kept
- std::vector<types::boundary_id> 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<unsigned int>
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// 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 <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/tria.h>
+
+#include "../tests.h"
+
+template <int dim>
+void
+print_boundary_ids(const dealii::Triangulation<dim> &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<dim> tria1;
+ Triangulation<dim> tria2;
+ Triangulation<dim> tria3;
+ Triangulation<dim> 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;
+}
--- /dev/null
+
+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::