]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Added Triangulation::create_cell_iterator(). 11365/head
authorMarc Fehling <mafehling.git@gmail.com>
Sun, 13 Dec 2020 02:41:12 +0000 (19:41 -0700)
committerMarc Fehling <mafehling.git@gmail.com>
Sun, 20 Dec 2020 21:08:57 +0000 (14:08 -0700)
12 files changed:
doc/news/changes/incompatibilities/20201212Fehling [new file with mode: 0644]
include/deal.II/grid/cell_id.h
include/deal.II/grid/tria.h
source/grid/cell_id.cc
source/grid/tria.cc
tests/fe/fe_subface_values_01.cc
tests/fullydistributed_grids/copy_distributed_tria_03.cc
tests/fullydistributed_grids/copy_serial_tria_08.cc
tests/fullydistributed_grids/copy_serial_tria_09.cc
tests/grid/cell_id_03.cc
tests/grid/cell_id_04.cc
tests/sharedtria/mg_dof_01.cc

diff --git a/doc/news/changes/incompatibilities/20201212Fehling b/doc/news/changes/incompatibilities/20201212Fehling
new file mode 100644 (file)
index 0000000..e90b30f
--- /dev/null
@@ -0,0 +1,4 @@
+Deprecated: The function CellId::to_cell() has been replaced by
+Triangulation::create_cell_iterator().
+<br>
+(Marc Fehling, 2020/12/12)
index 4b23bb371666683ab2016358204dfd8b3cc19da1..d1716bb6a4204fdfe3c0eb33b0a58da4b6082d0d 100644 (file)
@@ -148,9 +148,11 @@ public:
 
   /**
    * 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;
 
   /**
index bca79fbf5f33a2c69a2ed418cf22708eaa493c39..d98377b26ea5ad4c3f26801f5fcb5cfe88bb6f6b 100644 (file)
@@ -25,6 +25,7 @@
 #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>
@@ -2743,6 +2744,21 @@ public:
   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
    */
@@ -4029,8 +4045,6 @@ private:
 
   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
index 81ba05cdb629415ac89c7f3626f00b430ebb9963..2dbe9567e0e8b31bce160e52cf9f7864ea0ce903 100644 (file)
@@ -166,15 +166,10 @@ template <int dim, int spacedim>
 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"
 
index d9d2b7737ae7eeab3ec49fcc45df63fb639ca370..1feb29179cf11490c0a5c38399e2cf1f05fc4f6b 100644 (file)
@@ -11985,6 +11985,31 @@ Triangulation<dim, spacedim>::last_active() const
 
 
 
+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
index 195062eac059ed2c188a5cd9c49ae40c2d8491d9..fc9cc8f641a176534de62d56500cce2a6c139272 100644 (file)
@@ -74,12 +74,12 @@ run(unsigned int degree, unsigned int n_q_points)
   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;
index 30c067f4f8e219ef7cb9027935d113635c76922a..b420846261704d52a30a5e15ef7b8013d6ca55e7 100644 (file)
@@ -63,7 +63,7 @@ test(int n_refinements, const int n_subdivisions, MPI_Comm comm)
   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"));
index a83aa5ee74b7f0c1075d362e19ff0c1580b4b0e5..d0d68e5464b7318bd010ee6a6570b9d6cbda1637 100644 (file)
@@ -65,7 +65,7 @@ test(int n_refinements, MPI_Comm comm)
   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++)
index 118e2e3b7ea9f2d65752d31536545de09804a6e2..5bc44455ff2e5779c06ac8119396752f293c5010 100644 (file)
@@ -79,7 +79,7 @@ test(int                                              n_refinements,
   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"));
index 4d719abef022e13e9a723f3baa69249bc7b03223..8defc8c542d91bf8d73dbb50bdf0322174d92f4c 100644 (file)
@@ -49,7 +49,7 @@ check(TRIA &tr)
 
       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;
     }
index 55d79dc461a67fd11a687c74ee7ca7d4adae58cc..d66c42ccadeeb6a6cf13f7408f4b6bb934d8c601 100644 (file)
@@ -56,7 +56,8 @@ check(Triangulation<dim> &tr)
 
       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());
index 9918d477c2413678877d3307f8e2f4dfac4b698c..3fd2c10f23cb83cd362716a9eb579c339a705358 100644 (file)
@@ -75,7 +75,8 @@ compare_meshes(DoFHandler<dim> &shared_dof_handler,
             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(),

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.