// - always raise levels to match criterion, never lower them
// - exchange level indices on ghost cells
+ // Function that updates the level of neighbor to fulfill difference
+ // criterion, and returns whether it was changed.
+ const auto update_neighbor_level =
+ [&future_levels, max_difference](const active_cell_iterator &neighbor,
+ const level_type cell_level) -> bool {
+ // We only care about locally owned neighbors. If neighbor is a ghost cell,
+ // its future FE index will be updated on the owning process and
+ // communicated at the next loop iteration.
+ if (neighbor->is_locally_owned())
+ {
+ const level_type neighbor_level = static_cast<level_type>(
+ future_levels[neighbor->global_active_cell_index()]);
+
+ // ignore neighbors that are not part of the hierarchy
+ if (neighbor_level == invalid_level)
+ return false;
+
+ if ((cell_level - max_difference) > neighbor_level)
+ {
+ future_levels[neighbor->global_active_cell_index()] =
+ cell_level - max_difference;
+
+ return true;
+ }
+ }
+
+ return false;
+ };
+
bool levels_changed = false;
bool levels_changed_in_cycle;
do
for (unsigned int f = 0; f < cell->n_faces(); ++f)
if (cell->face(f)->at_boundary() == false)
{
- const auto neighbor = cell->neighbor(f);
-
- // We only care about locally owned neighbors. If neighbor is
- // a ghost cell, its future FE index will be updated on the
- // owning process and communicated at the next loop iteration.
- if (neighbor->is_locally_owned())
+ if (cell->face(f)->has_children())
{
- const level_type neighbor_level = static_cast<level_type>(
- future_levels[neighbor->global_active_cell_index()]);
-
- // ignore neighbors that are not part of the hierarchy
- if (neighbor_level == invalid_level)
- continue;
-
- if ((cell_level - max_difference) > neighbor_level)
+ for (unsigned int sf = 0;
+ sf < cell->face(f)->n_children();
+ ++sf)
{
- // update future level
- future_levels[neighbor->global_active_cell_index()] =
- cell_level - max_difference;
-
- levels_changed_in_cycle = true;
+ const auto neighbor =
+ cell->neighbor_child_on_subface(f, sf);
+ levels_changed_in_cycle |=
+ update_neighbor_level(neighbor, cell_level);
}
}
+ else
+ {
+ const auto neighbor = cell->neighbor(f);
+ levels_changed_in_cycle |=
+ update_neighbor_level(neighbor, cell_level);
+ }
}
}
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2021 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 restrictions on level differences imposed by
+// DoFHandler::prepare_coarsening_and_refinement()
+// on h-refined grids
+
+
+#include <deal.II/base/point.h>
+
+#include <deal.II/dofs/dof_handler.h>
+
+#include <deal.II/fe/fe_q.h>
+
+#include <deal.II/grid/tria.h>
+
+#include <deal.II/hp/fe_collection.h>
+
+#include "../tests.h"
+
+#include "../test_grids.h"
+
+
+unsigned int
+divide_and_ceil(unsigned int x, unsigned int y)
+{
+ return x / y + (x % y > 0);
+}
+
+
+template <int dim>
+void
+test(const unsigned int fes_size, const unsigned int max_difference)
+{
+ Assert(max_difference > 0, ExcInternalError());
+ Assert(fes_size > 1, ExcInternalError());
+
+ // setup FE collection
+ hp::FECollection<dim> fes;
+ while (fes.size() < fes_size)
+ fes.push_back(FE_Q<dim>(1));
+
+ const unsigned int contains_fe_index = 0;
+ const auto sequence = fes.get_hierarchy_sequence(contains_fe_index);
+
+ // set up line grid
+ // - refine leftmost column of cells consecutively
+ // - assign highest p-level to rightmost cell
+ //
+ // +++-+---+------+
+ // +++ | | |
+ // +++-+---+ |
+ // +++ | | |
+ // +++-+---+------+
+ //
+ // after prepare_coarsening_and_refinement(), each p-level will correspond to
+ // a unique h-level
+
+ Triangulation<dim> tria;
+ TestGrids::hyper_line(tria, 2);
+ const unsigned int n_refinements =
+ divide_and_ceil(sequence.size() - 1, max_difference);
+ for (unsigned int i = 0; i < n_refinements; ++i)
+ {
+ for (const auto &cell : tria.active_cell_iterators())
+ for (unsigned int v = 0; v < cell->n_vertices(); ++v)
+ if (cell->vertex(v)[0] == 0.)
+ cell->set_refine_flag();
+
+ tria.execute_coarsening_and_refinement();
+ }
+
+ DoFHandler<dim> dofh(tria);
+ for (const auto &cell : dofh.cell_iterators_on_level(0))
+ if (cell->is_active())
+ cell->set_active_fe_index(sequence.back());
+ dofh.distribute_dofs(fes);
+
+ const bool fe_indices_changed =
+ dofh.prepare_coarsening_and_refinement(max_difference, contains_fe_index);
+ tria.execute_coarsening_and_refinement();
+
+ (void)fe_indices_changed;
+ Assert(fe_indices_changed, ExcInternalError());
+
+#ifdef DEBUG
+ // check each cell's active FE by its h-level
+ for (unsigned int l = 0; l < tria.n_levels(); ++l)
+ for (const auto &cell : dofh.cell_iterators_on_level(l))
+ if (cell->is_active())
+ {
+ const unsigned int expected_level = std::max(
+ 0, static_cast<int>(sequence.size() - 1 - l * max_difference));
+ Assert(cell->active_fe_index() == sequence[expected_level],
+ ExcInternalError());
+ }
+#endif
+
+ deallog << "OK" << std::endl;
+}
+
+
+int
+main()
+{
+ initlog();
+
+ test<2>(5, 1);
+ test<2>(10, 2);
+ test<2>(15, 3);
+}
--- /dev/null
+
+DEAL::OK
+DEAL::OK
+DEAL::OK
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2021 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 restrictions on level differences imposed by
+// DoFHandler::prepare_coarsening_and_refinement()
+// on h-refined grids
+
+
+#include <deal.II/base/point.h>
+
+#include <deal.II/distributed/shared_tria.h>
+#include <deal.II/distributed/tria.h>
+#include <deal.II/distributed/tria_base.h>
+
+#include <deal.II/dofs/dof_handler.h>
+
+#include <deal.II/fe/fe_q.h>
+
+#include <deal.II/hp/fe_collection.h>
+
+#include "../tests.h"
+
+#include "../test_grids.h"
+
+
+unsigned int
+divide_and_ceil(unsigned int x, unsigned int y)
+{
+ return x / y + (x % y > 0);
+}
+
+
+template <int dim>
+void
+test(parallel::TriangulationBase<dim> &tria,
+ const unsigned int fes_size,
+ const unsigned int max_difference)
+{
+ Assert(tria.n_levels() == 0, ExcInternalError());
+ Assert(max_difference > 0, ExcInternalError());
+ Assert(fes_size > 1, ExcInternalError());
+
+ // setup FE collection
+ hp::FECollection<dim> fes;
+ while (fes.size() < fes_size)
+ fes.push_back(FE_Q<dim>(1));
+
+ const unsigned int contains_fe_index = 0;
+ const auto sequence = fes.get_hierarchy_sequence(contains_fe_index);
+
+ // set up line grid
+ // - refine leftmost column of cells consecutively
+ // - assign highest p-level to rightmost cell
+ //
+ // +++-+---+------+
+ // +++ | | |
+ // +++-+---+ |
+ // +++ | | |
+ // +++-+---+------+
+ //
+ // after prepare_coarsening_and_refinement(), each p-level will correspond to
+ // a unique column of cells and thus h-level
+
+ TestGrids::hyper_line(tria, 2);
+ const unsigned int n_refinements =
+ divide_and_ceil(sequence.size() - 1, max_difference);
+ for (unsigned int i = 0; i < n_refinements; ++i)
+ {
+ for (const auto &cell : tria.active_cell_iterators())
+ if (cell->is_locally_owned())
+ for (unsigned int v = 0; v < cell->n_vertices(); ++v)
+ if (cell->vertex(v)[0] == 0.)
+ cell->set_refine_flag();
+
+ tria.execute_coarsening_and_refinement();
+ }
+
+ DoFHandler<dim> dofh(tria);
+ for (const auto &cell : dofh.cell_iterators_on_level(0))
+ if (cell->is_active() && cell->is_locally_owned())
+ cell->set_active_fe_index(sequence.back());
+ dofh.distribute_dofs(fes);
+
+ const bool fe_indices_changed =
+ dofh.prepare_coarsening_and_refinement(max_difference, contains_fe_index);
+ tria.execute_coarsening_and_refinement();
+
+ (void)fe_indices_changed;
+ Assert(fe_indices_changed, ExcInternalError());
+
+#ifdef DEBUG
+ // check each cell's active FE by its h-level
+ for (unsigned int l = 0; l < tria.n_global_levels(); ++l)
+ for (const auto &cell : dofh.cell_iterators_on_level(l))
+ if (cell->is_active() && cell->is_locally_owned())
+ {
+ const unsigned int expected_level = std::max(
+ 0, static_cast<int>(sequence.size() - 1 - l * max_difference));
+ Assert(cell->active_fe_index() == sequence[expected_level],
+ ExcInternalError());
+ }
+#endif
+
+ deallog << "OK" << std::endl;
+}
+
+
+int
+main(int argc, char *argv[])
+{
+ Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);
+
+ if (Utilities::MPI::this_mpi_process(MPI_COMM_WORLD) == 0)
+ initlog();
+
+ constexpr const unsigned int dim = 2;
+
+ deallog << "parallel::shared::Triangulation" << std::endl;
+ {
+ parallel::shared::Triangulation<dim> tria(MPI_COMM_WORLD);
+
+ test<dim>(tria, 5, 1);
+ tria.clear();
+ test<dim>(tria, 10, 2);
+ tria.clear();
+ test<dim>(tria, 15, 3);
+ }
+
+ deallog << "parallel::shared::Triangulation with artificial cells"
+ << std::endl;
+ {
+ parallel::shared::Triangulation<dim> tria(MPI_COMM_WORLD,
+ Triangulation<dim>::none,
+ /*allow_artificial_cells=*/true);
+
+ test<dim>(tria, 5, 1);
+ tria.clear();
+ test<dim>(tria, 10, 2);
+ tria.clear();
+ test<dim>(tria, 15, 3);
+ }
+
+ deallog << "parallel::distributed::Triangulation" << std::endl;
+ {
+ parallel::distributed::Triangulation<dim> tria(MPI_COMM_WORLD);
+
+ test<dim>(tria, 5, 1);
+ tria.clear();
+ test<dim>(tria, 10, 2);
+ tria.clear();
+ test<dim>(tria, 15, 3);
+ }
+}
--- /dev/null
+
+DEAL::parallel::shared::Triangulation
+DEAL::OK
+DEAL::OK
+DEAL::OK
+DEAL::parallel::shared::Triangulation with artificial cells
+DEAL::OK
+DEAL::OK
+DEAL::OK
+DEAL::parallel::distributed::Triangulation
+DEAL::OK
+DEAL::OK
+DEAL::OK
--- /dev/null
+
+DEAL::parallel::shared::Triangulation
+DEAL::OK
+DEAL::OK
+DEAL::OK
+DEAL::parallel::shared::Triangulation with artificial cells
+DEAL::OK
+DEAL::OK
+DEAL::OK
+DEAL::parallel::distributed::Triangulation
+DEAL::OK
+DEAL::OK
+DEAL::OK