--- /dev/null
+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)
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.
*
}
+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>>
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.
*/
+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
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());
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));
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;
{
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);
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())
{
// 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());
{
// 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();
}
// 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());
{
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 =
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// 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();
+}
--- /dev/null
+
+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