]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Added cell->child_iterators(). 11000/head
authorMarc Fehling <mafehling.git@gmail.com>
Tue, 6 Oct 2020 01:31:06 +0000 (19:31 -0600)
committerMarc Fehling <mafehling.git@gmail.com>
Tue, 6 Oct 2020 01:31:06 +0000 (19:31 -0600)
13 files changed:
doc/news/changes/minor/20201003Fehling [new file with mode: 0644]
include/deal.II/dofs/dof_accessor.h
include/deal.II/dofs/dof_accessor.templates.h
include/deal.II/grid/tria_accessor.h
include/deal.II/grid/tria_accessor.templates.h
include/deal.II/numerics/cell_data_transfer.templates.h
source/distributed/error_predictor.cc
source/dofs/dof_handler.cc
source/hp/refinement.cc
source/numerics/solution_transfer.cc
source/particles/particle_handler.cc
tests/grid/accessor_03.cc [new file with mode: 0644]
tests/grid/accessor_03.output [new file with mode: 0644]

diff --git a/doc/news/changes/minor/20201003Fehling b/doc/news/changes/minor/20201003Fehling
new file mode 100644 (file)
index 0000000..42ca8ad
--- /dev/null
@@ -0,0 +1,5 @@
+New: Helper functions CellAccessor::child_iterators() and
+DoFCellAccessor::child_iterators() which return iterators to children of
+a cell via `cell->child_iterators()`.
+<br>
+(Marc Fehling, 2020/10/03)
index 997e5befff881a923fb38ecb6a41c51766050a23..c382b741818b1524e02e71787354b0fbc2c8ae49 100644 (file)
@@ -1470,6 +1470,15 @@ public:
   TriaIterator<DoFCellAccessor<dimension_, space_dimension_, level_dof_access>>
   child(const unsigned int i) const;
 
+  /**
+   * Return an array of iterators to all children of this cell.
+   */
+  boost::container::small_vector<
+    TriaIterator<
+      DoFCellAccessor<dimension_, space_dimension_, level_dof_access>>,
+    GeometryInfo<dimension_>::max_children_per_cell>
+  child_iterators() const;
+
   /**
    * Return an iterator to the @p ith face of this cell.
    *
index 0a31333faa8b7df51367569cb287d43355e21c03..3284b305fbe5318101188d7dd987312dd6ab632f 100644 (file)
@@ -2445,6 +2445,26 @@ DoFCellAccessor<dimension_, space_dimension_, level_dof_access>::child(
 }
 
 
+template <int dimension_, int space_dimension_, bool level_dof_access>
+inline boost::container::small_vector<
+  TriaIterator<DoFCellAccessor<dimension_, space_dimension_, level_dof_access>>,
+  GeometryInfo<dimension_>::max_children_per_cell>
+DoFCellAccessor<dimension_, space_dimension_, level_dof_access>::
+  child_iterators() const
+{
+  boost::container::small_vector<
+    TriaIterator<
+      DoFCellAccessor<dimension_, space_dimension_, level_dof_access>>,
+    GeometryInfo<dimension_>::max_children_per_cell>
+    child_iterators(this->n_children());
+
+  for (unsigned int i = 0; i < this->n_children(); ++i)
+    child_iterators[i] = this->child(i);
+
+  return child_iterators;
+}
+
+
 template <int dimension_, int space_dimension_, bool level_dof_access>
 inline TriaIterator<
   DoFCellAccessor<dimension_, space_dimension_, level_dof_access>>
index 3af0cabae73f1b8c6d1b5c8882c33d94baf86995..0cedfaf5826b8180d38d193eba9df78d2c7f372d 100644 (file)
@@ -2838,6 +2838,13 @@ public:
   TriaIterator<CellAccessor<dim, spacedim>>
   child(const unsigned int i) const;
 
+  /**
+   * Return an array of iterators to all children of this cell.
+   */
+  boost::container::small_vector<TriaIterator<CellAccessor<dim, spacedim>>,
+                                 GeometryInfo<dim>::max_children_per_cell>
+  child_iterators() const;
+
   /**
    * Return an iterator to the @p ith face of this cell.
    */
index e38575a64a80a031a4f6387c1917806e5aba13bd..4ab14d9b6abecb4d22500f89208feba1218f50b3 100644 (file)
@@ -3159,6 +3159,23 @@ CellAccessor<dim, spacedim>::child(const unsigned int i) const
 
 
 
+template <int dim, int spacedim>
+inline boost::container::small_vector<TriaIterator<CellAccessor<dim, spacedim>>,
+                                      GeometryInfo<dim>::max_children_per_cell>
+CellAccessor<dim, spacedim>::child_iterators() const
+{
+  boost::container::small_vector<TriaIterator<CellAccessor<dim, spacedim>>,
+                                 GeometryInfo<dim>::max_children_per_cell>
+    child_iterators(this->n_children());
+
+  for (unsigned int i = 0; i < this->n_children(); ++i)
+    child_iterators[i] = this->child(i);
+
+  return child_iterators;
+}
+
+
+
 template <int dim, int spacedim>
 inline TriaIterator<TriaAccessor<dim - 1, dim, spacedim>>
 CellAccessor<dim, spacedim>::face(const unsigned int i) const
index 91deef81734f85e2a35bfc12050487e2894ee7b0..14dde6d97ea2a631682fd73db3f926d3e3be506f 100644 (file)
@@ -105,11 +105,8 @@ CellDataTransfer<dim, spacedim, VectorType>::
               coarsened_cells_active_index.end())
             {
               std::set<unsigned int> indices_children;
-              for (unsigned int child_index = 0;
-                   child_index < parent->n_children();
-                   ++child_index)
+              for (const auto &sibling : parent->child_iterators())
                 {
-                  const auto &sibling = parent->child(child_index);
                   Assert(sibling->is_active() && sibling->coarsen_flag_set(),
                          typename dealii::Triangulation<
                            dim>::ExcInconsistentCoarseningFlags());
index d38da55a559e603afe4f113d1bc0a0f03465c716..ef12125b818aeb50d30da79ced9ac8f895d5788f 100644 (file)
@@ -258,12 +258,8 @@ namespace parallel
                 float sqrsum_of_predicted_errors = 0.;
                 float predicted_error            = 0.;
                 int   degree_difference          = 0;
-                for (unsigned int child_index = 0;
-                     child_index < cell->n_children();
-                     ++child_index)
+                for (const auto &child : cell->child_iterators())
                   {
-                    const auto child = cell->child(child_index);
-
                     predicted_error =
                       (**estimated_error_it)[child->active_cell_index()] /
                       (gamma_h * std::pow(.5, future_fe_degree));
@@ -336,10 +332,8 @@ namespace parallel
 
             case parallel::distributed::Triangulation<dim,
                                                       spacedim>::CELL_REFINE:
-              for (unsigned int child_index = 0;
-                   child_index < cell->n_children();
-                   ++child_index)
-                (**it_output)[cell->child(child_index)->active_cell_index()] =
+              for (const auto &child : cell->child_iterators())
+                (**it_output)[child->active_cell_index()] =
                   (*it_input) / std::sqrt(cell->n_children());
               break;
 
index 1335bd6f4cb1ef15ff321e74c73eaaaebfa9891d..22a0c26f5acb848c048948fdffd9209e5bf71c5b 100644 (file)
@@ -1951,11 +1951,8 @@ namespace internal
             {
               const auto &parent = refine.first;
 
-              for (unsigned int child_index = 0;
-                   child_index < parent->n_children();
-                   ++child_index)
+              for (const auto &child : parent->child_iterators())
                 {
-                  const auto &child = parent->child(child_index);
                   Assert(child->is_locally_owned() && child->is_active(),
                          ExcInternalError());
                   child->set_active_fe_index(refine.second);
index f287b9a9b16f9e05b3c4f43a81045a975f02bd90..3e80a37954317f0454e67ae08abb2374d5275d5d 100644 (file)
@@ -679,10 +679,8 @@ namespace hp
                 const unsigned int n_children = parent->n_children();
 
                 unsigned int h_flagged_children = 0, p_flagged_children = 0;
-                for (unsigned int child_index = 0; child_index < n_children;
-                     ++child_index)
+                for (const auto &child : parent->child_iterators())
                   {
-                    const auto &child = parent->child(child_index);
                     if (child->is_active())
                       {
                         if (child->is_locally_owned())
@@ -716,10 +714,8 @@ namespace hp
                   {
                     // Perform pure h coarsening and
                     // drop all p adaptation flags.
-                    for (unsigned int child_index = 0; child_index < n_children;
-                         ++child_index)
+                    for (const auto &child : parent->child_iterators())
                       {
-                        const auto &child = parent->child(child_index);
                         // h_flagged_children == n_children implies
                         // that all children are active
                         Assert(child->is_active(), ExcInternalError());
@@ -731,10 +727,8 @@ namespace hp
                   {
                     // Perform p adaptation on all children and
                     // drop all h coarsening flags.
-                    for (unsigned int child_index = 0; child_index < n_children;
-                         ++child_index)
+                    for (const auto &child : parent->child_iterators())
                       {
-                        const auto &child = parent->child(child_index);
                         if (child->is_active() && child->is_locally_owned())
                           child->clear_coarsen_flag();
                       }
index 9766483eebee2d0da6012ab486e9cf79bc8d84cb..67c25059ae9550be4be1bd0f701a828723a33d40 100644 (file)
@@ -357,10 +357,8 @@ SolutionTransfer<dim, VectorType, DoFHandlerType>::
           // case. we choose the 'least dominant fe' on all children from
           // the associated FECollection.
           std::set<unsigned int> fe_indices_children;
-          for (unsigned int child_index = 0; child_index < cell->n_children();
-               ++child_index)
+          for (const auto &child : cell->child_iterators())
             {
-              const auto &child = cell->child(child_index);
               Assert(child->is_active() && child->coarsen_flag_set(),
                      typename dealii::Triangulation<
                        dim>::ExcInconsistentCoarseningFlags());
index dda0e7b6c61007dc69e2d109d7dfc7e436237feb..a304691ffcead7a23e9528df3cd2622a5067713c 100644 (file)
@@ -1693,23 +1693,15 @@ namespace Particles
           {
             unsigned int n_particles = 0;
 
-            for (unsigned int child_index = 0;
-                 child_index < GeometryInfo<dim>::max_children_per_cell;
-                 ++child_index)
+            for (const auto &child : cell->child_iterators())
               {
-                const typename Triangulation<dim, spacedim>::cell_iterator
-                  child = cell->child(child_index);
                 n_particles += n_particles_in_cell(child);
               }
 
             stored_particles_on_cell.reserve(n_particles);
 
-            for (unsigned int child_index = 0;
-                 child_index < GeometryInfo<dim>::max_children_per_cell;
-                 ++child_index)
+            for (const auto &child : cell->child_iterators())
               {
-                const typename Triangulation<dim, spacedim>::cell_iterator
-                                         child       = cell->child(child_index);
                 const internal::LevelInd level_index = {child->level(),
                                                         child->index()};
                 const auto               particles_in_cell =
diff --git a/tests/grid/accessor_03.cc b/tests/grid/accessor_03.cc
new file mode 100644 (file)
index 0000000..2e779af
--- /dev/null
@@ -0,0 +1,74 @@
+// ---------------------------------------------------------------------
+//
+// 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.
+//
+// ---------------------------------------------------------------------
+
+
+
+// Verify functionality of CellAccessor::child_iterators() and
+// DoFCellAccessor:child_iterators()
+
+#include <deal.II/base/logstream.h>
+
+#include <deal.II/dofs/dof_accessor.h>
+#include <deal.II/dofs/dof_handler.h>
+
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/tria.h>
+#include <deal.II/grid/tria_accessor.h>
+
+#include "../tests.h"
+
+
+
+template <int dim>
+void
+test()
+{
+  Triangulation<dim> tria;
+  GridGenerator::hyper_cube(tria);
+  tria.refine_global(1);
+
+  {
+    const auto &parent = tria.begin(/*level=*/0);
+
+    for (const auto &child : parent->child_iterators())
+      deallog << child->id().to_string() << std::endl;
+  }
+
+  {
+    DoFHandler<dim> dofh(tria);
+
+    const auto &parent = dofh.begin(/*level=*/0);
+
+    for (const auto &child : parent->child_iterators())
+      deallog << child->id().to_string() << std::endl;
+  }
+}
+
+
+int
+main()
+{
+  initlog();
+
+  deallog.push("1d");
+  test<1>();
+  deallog.pop();
+  deallog.push("2d");
+  test<2>();
+  deallog.pop();
+  deallog.push("3d");
+  test<3>();
+  deallog.pop();
+}
diff --git a/tests/grid/accessor_03.output b/tests/grid/accessor_03.output
new file mode 100644 (file)
index 0000000..649733c
--- /dev/null
@@ -0,0 +1,29 @@
+
+DEAL:1d::0_1:0
+DEAL:1d::0_1:1
+DEAL:1d::0_1:0
+DEAL:1d::0_1:1
+DEAL:2d::0_1:0
+DEAL:2d::0_1:1
+DEAL:2d::0_1:2
+DEAL:2d::0_1:3
+DEAL:2d::0_1:0
+DEAL:2d::0_1:1
+DEAL:2d::0_1:2
+DEAL:2d::0_1:3
+DEAL:3d::0_1:0
+DEAL:3d::0_1:1
+DEAL:3d::0_1:2
+DEAL:3d::0_1:3
+DEAL:3d::0_1:4
+DEAL:3d::0_1:5
+DEAL:3d::0_1:6
+DEAL:3d::0_1:7
+DEAL:3d::0_1:0
+DEAL:3d::0_1:1
+DEAL:3d::0_1:2
+DEAL:3d::0_1:3
+DEAL:3d::0_1:4
+DEAL:3d::0_1:5
+DEAL:3d::0_1:6
+DEAL:3d::0_1:7

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.