]> https://gitweb.dealii.org/ - dealii.git/commitdiff
hide an implementational detail
authorMatthias Maier <tamiko@43-1.org>
Mon, 25 Nov 2019 22:41:02 +0000 (16:41 -0600)
committerMatthias Maier <tamiko@43-1.org>
Mon, 25 Nov 2019 23:33:02 +0000 (17:33 -0600)
The whole notion of a "Triangulation::active_cell_iterator" in DataOut and
DataOut_DoFData is unnecessary. This data type is in fact never used
except for one implementational detail.

Thus, remove it.

doc/news/changes/incompatibilities/20191125MatthiasMaier [new file with mode: 0644]
include/deal.II/numerics/data_out.h
include/deal.II/numerics/data_out_dof_data.h
include/deal.II/numerics/data_out_dof_data.templates.h
source/numerics/data_out.cc

diff --git a/doc/news/changes/incompatibilities/20191125MatthiasMaier b/doc/news/changes/incompatibilities/20191125MatthiasMaier
new file mode 100644 (file)
index 0000000..eb0152e
--- /dev/null
@@ -0,0 +1,5 @@
+Removed: The typedefs DataOut::active_cell_iterator and
+DataOut_DoFData::active_cell_iterator have been removed. These typedefs are
+an implementational detail that should have never exposed in the interface.
+<br>
+(Matthias Maier, 2019/11/25)
index e0051930d491649d8eed2cb9be88cd42de90a672..0d57f5c5a03b40b3a134bb1001f1124e284a42bd 100644 (file)
@@ -170,10 +170,6 @@ public:
     typename DataOut_DoFData<DoFHandlerType,
                              DoFHandlerType::dimension,
                              DoFHandlerType::space_dimension>::cell_iterator;
-  using active_cell_iterator = typename DataOut_DoFData<
-    DoFHandlerType,
-    DoFHandlerType::dimension,
-    DoFHandlerType::space_dimension>::active_cell_iterator;
 
   /**
    * Enumeration describing the part of the domain in which cells
index 7d6622dac45c87c6e7f1fe7f775867ba7004dde3..0c5669ee36f75913446ce1db53f162a06a2b36ce 100644 (file)
@@ -604,9 +604,6 @@ public:
   using cell_iterator =
     typename Triangulation<DoFHandlerType::dimension,
                            DoFHandlerType::space_dimension>::cell_iterator;
-  using active_cell_iterator = typename Triangulation<
-    DoFHandlerType::dimension,
-    DoFHandlerType::space_dimension>::active_cell_iterator;
 
 public:
   /**
index 3675e60d46f0d3d1ad27ccee903a23a0d3dec0a9..c982eede9c002da9d34edd962521f08c8a40ac4a 100644 (file)
@@ -228,7 +228,11 @@ namespace internal
           bool duplicate = false;
           for (unsigned int j = 0; j < dataset; ++j)
             if (finite_elements[dataset].get() == finite_elements[j].get())
-              duplicate = true;
+              {
+                duplicate = true;
+                break;
+              }
+
           if (duplicate == false)
             {
               typename DoFHandlerType::active_cell_iterator dh_cell(
index 96f77701fc94ccfdf22cf14eba60ad8b5ddf754e..eac86627f28c57003fe0e3f50c1445127724534f 100644 (file)
@@ -910,16 +910,16 @@ DataOut<dim, DoFHandlerType>::build_patches(
     // data from (cell data vectors do not have the length distance computed by
     // first_locally_owned_cell/next_locally_owned_cell because this might skip
     // some values (FilteredIterator).
-    active_cell_iterator active_cell  = this->triangulation->begin_active();
-    unsigned int         active_index = 0;
-    cell_iterator        cell         = first_locally_owned_cell();
+    auto          active_cell  = this->triangulation->begin_active();
+    unsigned int  active_index = 0;
+    cell_iterator cell         = first_locally_owned_cell();
     for (; cell != this->triangulation->end();
          cell = next_locally_owned_cell(cell))
       {
         // move forward until active_cell points at the cell (cell) we are
         // looking at to compute the current active_index
         while (active_cell != this->triangulation->end() && cell->active() &&
-               active_cell_iterator(cell) != active_cell)
+               decltype(active_cell)(cell) != active_cell)
           {
             ++active_cell;
             ++active_index;

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.