std::pair<std::pair<cell_iterator, unsigned int>, std::bitset<3>>> &
get_periodic_face_map() const;
+ /**
+ * Return vector filled with the used reference-cell types of this
+ * triangulation.
+ */
+ const std::vector<ReferenceCell::Type> &
+ get_reference_cell_types() const;
+
+ /**
+ * Indicate if the triangulation only consists of hypercube-like cells, i.e.,
+ * lines, quadrilaterals, or hexahedrons.
+ */
+ bool
+ all_reference_cell_types_are_hyper_cube() const;
+
#ifdef DOXYGEN
/**
* Write and read the data of this object from a stream for the purpose
*/
MeshSmoothing smooth_grid;
+ /**
+ * Vector caching all reference-cell types of the given triangulation
+ * (also in the distributed case).
+ */
+ std::vector<ReferenceCell::Type> reference_cell_types;
+
/**
* Write a bool vector to the given stream, writing a pre- and a postfix
* magic number. The vector is written in an almost binary format, i.e. the
void
update_periodic_face_map();
+ /**
+ * Update the internal reference_cell_types vector.
+ */
+ virtual void
+ update_reference_cell_types();
+
private:
/**
#endif
+ template <int dim, int spacedim>
+ void
+ TriangulationBase<dim, spacedim>::update_reference_cell_types()
+ {
+ // run algorithm for locally-owned cells
+ dealii::Triangulation<dim, spacedim>::update_reference_cell_types();
+
+ // translate ReferenceCell::Type to unsigned int (needed by
+ // Utilities::MPI::compute_set_union)
+ std::vector<unsigned int> reference_cell_types_ui;
+
+ for (const auto &i : this->reference_cell_types)
+ reference_cell_types_ui.push_back(static_cast<unsigned int>(i));
+
+ // create union
+ reference_cell_types_ui =
+ Utilities::MPI::compute_set_union(reference_cell_types_ui,
+ this->mpi_communicator);
+
+ // transform back and store result
+ this->reference_cell_types.clear();
+ for (const auto &i : reference_cell_types_ui)
+ this->reference_cell_types.push_back(static_cast<ReferenceCell::Type>(i));
+ }
+
+
template <int dim, int spacedim>
types::subdomain_id
TriangulationBase<dim, spacedim>::locally_owned_subdomain() const
Triangulation<dim, spacedim> &&tria) noexcept
: Subscriptor(std::move(tria))
, smooth_grid(tria.smooth_grid)
+ , reference_cell_types(std::move(tria.reference_cell_types))
, periodic_face_pairs_level_0(std::move(tria.periodic_face_pairs_level_0))
, periodic_face_map(std::move(tria.periodic_face_map))
, levels(std::move(tria.levels))
Subscriptor::operator=(std::move(tria));
smooth_grid = tria.smooth_grid;
+ reference_cell_types = std::move(tria.reference_cell_types);
periodic_face_pairs_level_0 = std::move(tria.periodic_face_pairs_level_0);
periodic_face_map = std::move(tria.periodic_face_map);
levels = std::move(tria.levels);
clear_despite_subscriptions();
periodic_face_pairs_level_0.clear();
periodic_face_map.clear();
+ reference_cell_types.clear();
}
vertices_used = other_tria.vertices_used;
anisotropic_refinement = other_tria.anisotropic_refinement;
smooth_grid = other_tria.smooth_grid;
+ reference_cell_types = other_tria.reference_cell_types;
if (dim > 1)
faces = std::make_unique<internal::TriangulationImplementation::TriaFaces>(
// because sometimes other objects are already attached to it:
try
{
-#ifndef DEAL_II_WITH_SIMPLEX_SUPPORT
- AssertThrow(
- std::any_of(cells.begin(),
- cells.end(),
- [](const auto &cell) {
- return cell.vertices.size() !=
- GeometryInfo<dim>::vertices_per_cell;
- }) == false,
- ExcMessage(
- "A cell with invalid number of vertices has been provided."));
-#endif
-
internal::TriangulationImplementation::Implementation::
create_triangulation(v, cells, subcelldata, *this);
+
+ this->update_reference_cell_types();
+
+#ifndef DEAL_II_WITH_SIMPLEX_SUPPORT
+ Assert(this->all_reference_cell_types_are_hyper_cube(),
+ ExcMessage(
+ "A cell with invalid number of vertices has been provided."));
+#endif
}
catch (...)
{
+template <int dim, int spacedim>
+void
+Triangulation<dim, spacedim>::update_reference_cell_types()
+{
+ std::set<ReferenceCell::Type> reference_cell_types_set;
+ for (auto cell : active_cell_iterators())
+ if (cell->is_locally_owned())
+ reference_cell_types_set.insert(cell->reference_cell_type());
+
+ std::vector<ReferenceCell::Type> reference_cell_types(
+ reference_cell_types_set.begin(), reference_cell_types_set.end());
+
+ this->reference_cell_types = reference_cell_types;
+}
+
+
+
+template <int dim, int spacedim>
+const std::vector<ReferenceCell::Type> &
+Triangulation<dim, spacedim>::get_reference_cell_types() const
+{
+ return this->reference_cell_types;
+}
+
+
+
+template <int dim, int spacedim>
+bool
+Triangulation<dim, spacedim>::all_reference_cell_types_are_hyper_cube() const
+{
+ return (this->reference_cell_types.size() == 0) ||
+ (this->reference_cell_types.size() == 1 &&
+ this->reference_cell_types[0] == ReferenceCell::get_hypercube(dim));
+}
+
+
+
template <int dim, int spacedim>
void
Triangulation<dim, spacedim>::clear_despite_subscriptions()
}
+
template <int dim, int spacedim>
typename Triangulation<dim, spacedim>::DistortedCellList
Triangulation<dim, spacedim>::execute_refinement()