virtual std::map<unsigned int, std::set<dealii::types::subdomain_id>>
compute_vertices_with_ghost_neighbors() const;
+ /**
+ * @copydoc dealii::Triangulation::get_boundary_ids()
+ *
+ * @note This function involves a global communication gathering all current
+ * IDs from all processes.
+ */
+ virtual std::vector<types::boundary_id>
+ get_boundary_ids() const override;
+
+ /**
+ * @copydoc dealii::Triangulation::get_manifold_ids()
+ *
+ * @note This function involves a global communication gathering all current
+ * IDs from all processes.
+ */
+ virtual std::vector<types::manifold_id>
+ get_manifold_ids() const override;
+
protected:
/**
* MPI communicator to be used for the triangulation. We create a unique
/**
* Return a vector containing all boundary indicators assigned to boundary
- * faces of this Triangulation object. Note, that each boundary indicator is
- * reported only once. The size of the return vector will represent the
- * number of different indicators (which is greater or equal one).
+ * faces of active cells of this Triangulation object. Note, that each
+ * boundary indicator is reported only once. The size of the return vector
+ * will represent the number of different indicators (which is greater or
+ * equal one).
*
* @ingroup boundary
*
* @see
* @ref GlossBoundaryIndicator "Glossary entry on boundary indicators"
*/
- std::vector<types::boundary_id>
+ virtual std::vector<types::boundary_id>
get_boundary_ids() const;
/**
* Return a vector containing all manifold indicators assigned to the
- * objects of this Triangulation. Note, that each manifold indicator is
- * reported only once. The size of the return vector will represent the
- * number of different indicators (which is greater or equal one).
+ * objects of the active cells of this Triangulation. Note, that each
+ * manifold indicator is reported only once. The size of the return vector
+ * will represent the number of different indicators (which is greater or
+ * equal one).
*
* @ingroup manifold
*
* @see
* @ref GlossManifoldIndicator "Glossary entry on manifold indicators"
*/
- std::vector<types::manifold_id>
+ virtual std::vector<types::manifold_id>
get_manifold_ids() const;
/**
+ template <int dim, int spacedim>
+ std::vector<types::boundary_id>
+ TriangulationBase<dim, spacedim>::get_boundary_ids() const
+ {
+ return Utilities::MPI::compute_set_union(
+ dealii::Triangulation<dim, spacedim>::get_boundary_ids(),
+ this->mpi_communicator);
+ }
+
+
+
+ template <int dim, int spacedim>
+ std::vector<types::manifold_id>
+ TriangulationBase<dim, spacedim>::get_manifold_ids() const
+ {
+ return Utilities::MPI::compute_set_union(
+ dealii::Triangulation<dim, spacedim>::get_manifold_ids(),
+ this->mpi_communicator);
+ }
+
+
+
template <int dim, int spacedim>
DistributedTriangulationBase<dim, spacedim>::DistributedTriangulationBase(
MPI_Comm mpi_communicator,
else
{
std::set<types::boundary_id> b_ids;
- active_cell_iterator cell = begin_active();
- for (; cell != end(); ++cell)
- for (const unsigned int face : GeometryInfo<dim>::face_indices())
- if (cell->at_boundary(face))
- b_ids.insert(cell->face(face)->boundary_id());
+ for (auto cell : active_cell_iterators())
+ if (cell->is_locally_owned())
+ for (const unsigned int face : GeometryInfo<dim>::face_indices())
+ if (cell->at_boundary(face))
+ b_ids.insert(cell->face(face)->boundary_id());
std::vector<types::boundary_id> boundary_ids(b_ids.begin(), b_ids.end());
return boundary_ids;
}
Triangulation<dim, spacedim>::get_manifold_ids() const
{
std::set<types::manifold_id> m_ids;
- active_cell_iterator cell = begin_active();
- for (; cell != end(); ++cell)
- {
- m_ids.insert(cell->manifold_id());
- if (dim > 1)
- for (const unsigned int face : GeometryInfo<dim>::face_indices())
- if (cell->at_boundary(face))
- m_ids.insert(cell->face(face)->manifold_id());
- }
+ for (auto cell : active_cell_iterators())
+ if (cell->is_locally_owned())
+ {
+ m_ids.insert(cell->manifold_id());
+ if (dim > 1)
+ for (const unsigned int face : GeometryInfo<dim>::face_indices())
+ if (cell->at_boundary(face))
+ m_ids.insert(cell->face(face)->manifold_id());
+ }
std::vector<types::manifold_id> manifold_indicators(m_ids.begin(),
m_ids.end());
return manifold_indicators;
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 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.
+//
+// ---------------------------------------------------------------------
+
+
+// create a tria mesh and copy it
+
+#include <deal.II/distributed/tria.h>
+
+#include <deal.II/grid/grid_generator.h>
+
+#include "../tests.h"
+
+template <int dim>
+void
+test()
+{
+ parallel::distributed::Triangulation<dim> tria(MPI_COMM_WORLD);
+ GridGenerator::subdivided_hyper_cube(tria, 4);
+
+ for (auto cell : tria.active_cell_iterators())
+ if (cell->is_locally_owned())
+ for (const unsigned int face : GeometryInfo<dim>::face_indices())
+ {
+ if (std::fabs(cell->face(face)->center()(0) - 0.0) < 1e-12)
+ cell->face(face)->set_all_boundary_ids(1);
+ if (std::fabs(cell->face(face)->center()(0) - 1.0) < 1e-12)
+ cell->face(face)->set_all_boundary_ids(2);
+ if (std::fabs(cell->face(face)->center()(1) - 0.0) < 1e-12)
+ cell->face(face)->set_all_boundary_ids(3);
+ if (std::fabs(cell->face(face)->center()(1) - 1.0) < 1e-12)
+ cell->face(face)->set_all_boundary_ids(4);
+ }
+
+ for (auto i : tria.get_boundary_ids())
+ deallog << i << " ";
+ deallog << std::endl;
+}
+
+int
+main(int argc, char *argv[])
+{
+ Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);
+ MPILogInitAll all;
+
+ test<2>();
+}
--- /dev/null
+
+DEAL:0::1 2 3 4
+
+DEAL:1::1 2 3 4
+