--- /dev/null
+Deprecated: The function CellId::to_cell() has been replaced by
+Triangulation::create_cell_iterator().
+<br>
+(Marc Fehling, 2020/12/12)
/**
* Return a cell_iterator to the cell represented by this CellId.
+ *
+ * @deprecated Use Triangulation::create_cell_iterator() instead.
*/
template <int dim, int spacedim>
- typename Triangulation<dim, spacedim>::cell_iterator
+ DEAL_II_DEPRECATED typename Triangulation<dim, spacedim>::cell_iterator
to_cell(const Triangulation<dim, spacedim> &tria) const;
/**
#include <deal.II/base/smartpointer.h>
#include <deal.II/base/subscriptor.h>
+#include <deal.II/grid/cell_id.h>
#include <deal.II/grid/tria_description.h>
#include <deal.II/grid/tria_iterator_selector.h>
#include <deal.II/grid/tria_levels.h>
active_cell_iterator
last_active() const;
+ /**
+ * Return an iterator to a cell of this Triangulation object constructed from
+ * an independent CellId object.
+ *
+ * If the given argument corresponds to a valid cell in this triangulation,
+ * this operation will always succeed for sequential triangulations where the
+ * current processor stores all cells that are part of the triangulation. On
+ * the other hand, if this is a parallel triangulation, then the current
+ * processor may not actually know about this cell. In this case, this
+ * operation will succeed for locally relevant cells, but may not for
+ * artificial cells that are less refined on the current processor.
+ */
+ cell_iterator
+ create_cell_iterator(const CellId &cell_id) const;
+
/**
* @name Cell iterator functions returning ranges of iterators
*/
friend class dealii::internal::TriangulationImplementation::TriaObjects;
- friend class CellId;
-
// explicitly check for sensible template arguments, but not on windows
// because MSVC creates bogus warnings during normal compilation
#ifndef DEAL_II_MSVC
typename Triangulation<dim, spacedim>::cell_iterator
CellId::to_cell(const Triangulation<dim, spacedim> &tria) const
{
- typename Triangulation<dim, spacedim>::cell_iterator cell(
- &tria, 0, tria.coarse_cell_id_to_coarse_cell_index(coarse_cell_id));
-
- for (unsigned int i = 0; i < n_child_indices; ++i)
- cell = cell->child(static_cast<unsigned int>(child_indices[i]));
-
- return cell;
+ return tria.create_cell_iterator(*this);
}
+
// explicit instantiations
#include "cell_id.inst"
+template <int dim, int spacedim>
+typename Triangulation<dim, spacedim>::cell_iterator
+Triangulation<dim, spacedim>::create_cell_iterator(const CellId &cell_id) const
+{
+ cell_iterator cell(
+ this, 0, coarse_cell_id_to_coarse_cell_index(cell_id.get_coarse_cell_id()));
+
+ for (const auto &child_index : cell_id.get_child_indices())
+ {
+ Assert(
+ cell->has_children(),
+ ExcMessage(
+ "CellId is invalid for this triangulation.\n"
+ "Either the provided CellId does not correspond to a cell in this "
+ "triangulation object, or, in case you are using a parallel "
+ "triangulation, may correspond to an artificial cell that is less "
+ "refined on this processor."));
+ cell = cell->child(static_cast<unsigned int>(child_index));
+ }
+
+ return cell;
+}
+
+
+
template <int dim, int spacedim>
typename Triangulation<dim, spacedim>::cell_iterator
Triangulation<dim, spacedim>::end() const
UpdateFlags flags = update_quadrature_points; // dummy
FESubfaceValues<dim> sub_face(fe, quad, flags);
- sub_face.reinit(CellId(0, {1}).to_cell(tria), 0, 0);
+ sub_face.reinit(tria.create_cell_iterator(CellId(0, {1})), 0, 0);
for (auto p : sub_face.get_quadrature_points())
deallog << p << std::endl;
deallog << std::endl;
- sub_face.reinit(CellId(0, {1}).to_cell(tria), 1, 0);
+ sub_face.reinit(tria.create_cell_iterator(CellId(0, {1})), 1, 0);
for (auto p : sub_face.get_quadrature_points())
deallog << p << std::endl;
deallog << std::endl;
for (auto &cell : tria_pft.active_cell_iterators())
{
CellId id = cell->id();
- auto cell_base = id.to_cell(tria_pdt);
+ auto cell_base = tria_pdt.create_cell_iterator(id);
for (unsigned int d = 0; d < dim; d++)
Assert(std::abs(cell->center()[d] - cell_base->center()[d]) < 1e-9,
ExcMessage("Cells do not match"));
for (auto &cell : tria_pft.active_cell_iterators())
{
CellId id = cell->id();
- auto cell_base = id.to_cell(basetria);
+ auto cell_base = basetria.create_cell_iterator(id);
// Assert(cell->center() == cell_base->center(),
// ExcMessage("Cells do not match"));
for (unsigned int d = 0; d < dim; d++)
for (auto &cell : tria_pft.active_cell_iterators())
{
CellId id = cell->id();
- auto cell_base = id.to_cell(basetria);
+ auto cell_base = basetria.create_cell_iterator(id);
for (unsigned int d = 0; d < dim; d++)
Assert(std::abs(cell->center()[d] - cell_base->center()[d]) < 1e-9,
ExcMessage("Cells do not match"));
const CellId cid = cell->id();
- typename TRIA::cell_iterator cell2 = cid.to_cell(tr);
+ typename TRIA::cell_iterator cell2 = tr.create_cell_iterator(cid);
deallog << cell2->level() << " " << cell2->index() << std::endl;
}
const CellId cid2(cids);
- typename Triangulation<dim>::cell_iterator cell2 = cid2.to_cell(tr);
+ typename Triangulation<dim>::cell_iterator cell2 =
+ tr.create_cell_iterator(cid2);
Assert(cid2 == cid, ExcInternalError());
Assert(cell2 == cell, ExcInternalError());
continue;
typename Triangulation<dim>::cell_iterator tria_shared_cell =
- cell->id().to_cell(shared_dof_handler.get_triangulation());
+ shared_dof_handler.get_triangulation().create_cell_iterator(
+ cell->id());
typename DoFHandler<dim>::cell_iterator dof_shared_cell(
&shared_dof_handler.get_triangulation(),
tria_shared_cell->level(),