Description<dim, spacedim>
create_description_from_triangulation(
const Triangulation<dim, spacedim> & tria,
- const LinearAlgebra::distributed::Vector<double> &partition)
+ const LinearAlgebra::distributed::Vector<double> &partition,
+ const TriangulationDescription::Settings settings)
+ {
+ const bool construct_multigrid =
+ (partition.size() > 0) &&
+ (settings &
+ TriangulationDescription::Settings::construct_multigrid_hierarchy);
+
+ Assert(
+ construct_multigrid == false ||
+ (tria.get_mesh_smoothing() &
+ Triangulation<dim, spacedim>::limit_level_difference_at_vertices),
+ ExcMessage(
+ "Source triangulation has to be set up with "
+ "limit_level_difference_at_vertices if the construction of the "
+ "multigrid hierarchy is requested!"));
+
+ std::vector<LinearAlgebra::distributed::Vector<double>> partitions_mg;
+
+ if (construct_multigrid) // perform first child policy
+ {
+ const auto tria_parallel =
+ dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(
+ &tria);
+
+ Assert(tria_parallel, ExcInternalError());
+
+ partition.update_ghost_values();
+
+ partitions_mg.resize(tria.n_global_levels());
+
+ for (unsigned int l = 0; l < tria.n_global_levels(); ++l)
+ partitions_mg[l].reinit(
+ tria_parallel->global_level_cell_index_partitioner(l).lock());
+
+ for (int level = tria.n_global_levels() - 1; level >= 0; --level)
+ {
+ for (const auto &cell : tria.cell_iterators_on_level(level))
+ {
+ if (cell->is_locally_owned_on_level() == false)
+ continue;
+
+ if (cell->is_active())
+ partitions_mg[level][cell->global_level_cell_index()] =
+ partition[cell->global_active_cell_index()];
+ else
+ partitions_mg[level][cell->global_level_cell_index()] =
+ partitions_mg[level + 1]
+ [cell->child(0)->global_level_cell_index()];
+ }
+
+ partitions_mg[level].update_ghost_values();
+ }
+ }
+
+ return create_description_from_triangulation(tria,
+ partition,
+ partitions_mg,
+ settings);
+ }
+
+
+
+ template <int dim, int spacedim>
+ Description<dim, spacedim>
+ create_description_from_triangulation(
+ const Triangulation<dim, spacedim> & tria,
+ const LinearAlgebra::distributed::Vector<double> &partition,
+ const std::vector<LinearAlgebra::distributed::Vector<double>>
+ & partitions_mg,
+ const TriangulationDescription::Settings settings_in)
{
#ifdef DEAL_II_WITH_MPI
if (tria.get_communicator() == MPI_COMM_NULL)
if (partition.size() == 0)
{
+ AssertDimension(partitions_mg.size(), 0);
return create_description_from_triangulation(tria,
- tria.get_communicator());
+ tria.get_communicator(),
+ settings_in);
}
partition.update_ghost_values();
+ for (const auto &partition : partitions_mg)
+ partition.update_ghost_values();
+
// 1) determine processes owning locally owned cells
const std::vector<unsigned int> relevant_processes = [&]() {
std::set<unsigned int> relevant_processes;
relevant_processes.insert(
static_cast<unsigned int>(partition.local_element(i)));
+ for (const auto &partition : partitions_mg)
+ for (unsigned int i = 0; i < partition.local_size(); ++i)
+ relevant_processes.insert(
+ static_cast<unsigned int>(partition.local_element(i)));
+
return std::vector<unsigned int>(relevant_processes.begin(),
relevant_processes.end());
}();
- TriangulationDescription::Settings settings =
- TriangulationDescription::Settings::default_setting;
+ const bool construct_multigrid = partitions_mg.size() > 0;
+
+ TriangulationDescription::Settings settings = settings_in;
+
+ if (construct_multigrid)
+ settings = static_cast<TriangulationDescription::Settings>(
+ settings |
+ TriangulationDescription::Settings::construct_multigrid_hierarchy);
const auto subdomain_id_function = [&partition](const auto &cell) {
if ((cell->is_active() && (cell->is_artificial() == false)))
return numbers::artificial_subdomain_id;
};
- const auto level_subdomain_id_function = [](const auto & /*cell*/) {
- return numbers::artificial_subdomain_id;
- };
+ const auto level_subdomain_id_function =
+ [&construct_multigrid, &partitions_mg](const auto &cell) {
+ if (construct_multigrid && (cell->is_artificial_on_level() == false))
+ return static_cast<unsigned int>(
+ partitions_mg[cell->level()][cell->global_level_cell_index()]);
+ else
+ return numbers::artificial_subdomain_id;
+ };
CreateDescriptionFromTriangulationHelper<dim, spacedim> helper(
tria,
--- /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.
+//
+// ---------------------------------------------------------------------
+
+
+// Test
+// TriangulationDescription::Utilities::create_description_from_triangulation()
+// so that it also works for p:d:T set up on a subcommunicator.
+
+#include <deal.II/base/mpi_consensus_algorithms.h>
+
+#include <deal.II/distributed/fully_distributed_tria.h>
+#include <deal.II/distributed/repartitioning_policy_tools.h>
+#include <deal.II/distributed/tria.h>
+
+#include <deal.II/fe/fe_q.h>
+
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/grid_out.h>
+#include <deal.II/grid/grid_tools.h>
+#include <deal.II/grid/tria_description.h>
+
+#include <deal.II/lac/la_parallel_vector.h>
+
+#include <deal.II/numerics/data_out.h>
+
+#include "./tests.h"
+
+using namespace dealii;
+
+namespace dealii
+{
+ namespace parallel
+ {
+ template <int dim, int spacedim = dim>
+ std::function<
+ unsigned int(const typename Triangulation<dim, spacedim>::cell_iterator &,
+ const typename Triangulation<dim, spacedim>::CellStatus)>
+ hanging_nodes_weighting(const double weight)
+ {
+ return [weight](const auto &cell, const auto &) -> unsigned int {
+ if (cell->is_locally_owned() == false)
+ return 0;
+
+ bool flag = false;
+
+ for (const auto f : cell->face_indices())
+ if (!cell->at_boundary(f) &&
+ (cell->neighbor(f)->has_children() ||
+ cell->level() != cell->neighbor(f)->level()))
+ flag = true;
+
+ if (flag)
+ return 10000 * weight;
+ else
+ return 10000;
+ };
+ }
+ } // namespace parallel
+} // namespace dealii
+
+
+template <int dim, int spacedim>
+void
+create_triangulation(Triangulation<dim, spacedim> &tria)
+{
+ const unsigned int n_ref_global = 3;
+ const unsigned int n_ref_local = 3;
+ GridGenerator::hyper_cube(tria, -1.0, +1.0);
+ tria.refine_global(n_ref_global);
+
+ for (unsigned int i = 0; i < n_ref_local; ++i)
+ {
+ for (auto cell : tria.active_cell_iterators())
+ if (cell->is_locally_owned())
+ {
+ bool flag = true;
+ for (unsigned int d = 0; d < dim; ++d)
+ if (cell->center()[d] > 0.0)
+ flag = false;
+ if (flag)
+ cell->set_refine_flag();
+ }
+ tria.execute_coarsening_and_refinement();
+ }
+}
+
+
+template <int dim>
+void
+test(const MPI_Comm comm)
+{
+ parallel::distributed::Triangulation<dim> tria(
+ comm,
+ Triangulation<dim>::none,
+ parallel::distributed::Triangulation<dim>::construct_multigrid_hierarchy);
+ create_triangulation(tria);
+
+ const RepartitioningPolicyTools::FirstChildPolicy<dim> policy(tria);
+
+ // const RepartitioningPolicyTools::CellWeightPolicy<dim>
+ // policy(parallel::hanging_nodes_weighting<dim>(2.0));
+
+ const auto construction_data =
+ TriangulationDescription::Utilities::create_description_from_triangulation(
+ tria,
+ policy.partition(tria),
+ TriangulationDescription::Settings::construct_multigrid_hierarchy);
+
+ parallel::fullydistributed::Triangulation<dim> tria_pft(comm);
+ tria_pft.create_triangulation(construction_data);
+
+ if (false)
+ {
+ GridOut grid_out;
+ grid_out.write_mesh_per_processor_as_vtu(tria_pft, "test", true, true);
+ }
+
+ FE_Q<dim> fe(2);
+
+ DoFHandler<dim> dof_handler(tria);
+ dof_handler.distribute_dofs(fe);
+ dof_handler.distribute_mg_dofs();
+
+ DoFHandler<dim> dof_handler_pft(tria_pft);
+ dof_handler_pft.distribute_dofs(fe);
+ dof_handler_pft.distribute_mg_dofs();
+
+ // print statistics
+ print_statistics(tria, true);
+ print_statistics(tria_pft, true);
+ print_statistics(dof_handler, true);
+ print_statistics(dof_handler_pft, true);
+}
+
+
+
+int
+main(int argc, char **argv)
+{
+ Utilities::MPI::MPI_InitFinalize mpi(argc, argv, 1);
+ MPILogInitAll all;
+
+ MPI_Comm comm = MPI_COMM_WORLD;
+
+ deallog.push("all");
+ test<2>(comm);
+ deallog.pop();
+}
--- /dev/null
+
+DEAL:0:all::n_levels: 7
+DEAL:0:all::n_cells: 577
+DEAL:0:all::n_active_cells: 433
+DEAL:0:all::n_cells on level=0: 1
+DEAL:0:all::n_cells on level=1: 4
+DEAL:0:all::n_cells on level=2: 16
+DEAL:0:all::n_cells on level=3: 24
+DEAL:0:all::n_cells on level=4: 44
+DEAL:0:all::n_cells on level=5: 116
+DEAL:0:all::n_cells on level=6: 372
+DEAL:0:all::n_active_cells on level=0: 0
+DEAL:0:all::n_active_cells on level=1: 0
+DEAL:0:all::n_active_cells on level=2: 10
+DEAL:0:all::n_active_cells on level=3: 13
+DEAL:0:all::n_active_cells on level=4: 15
+DEAL:0:all::n_active_cells on level=5: 23
+DEAL:0:all::n_active_cells on level=6: 372
+DEAL:0:all::
+DEAL:0:all::n_levels: 7
+DEAL:0:all::n_cells: 577
+DEAL:0:all::n_active_cells: 433
+DEAL:0:all::n_cells on level=0: 1
+DEAL:0:all::n_cells on level=1: 4
+DEAL:0:all::n_cells on level=2: 16
+DEAL:0:all::n_cells on level=3: 24
+DEAL:0:all::n_cells on level=4: 44
+DEAL:0:all::n_cells on level=5: 116
+DEAL:0:all::n_cells on level=6: 372
+DEAL:0:all::n_active_cells on level=0: 0
+DEAL:0:all::n_active_cells on level=1: 0
+DEAL:0:all::n_active_cells on level=2: 10
+DEAL:0:all::n_active_cells on level=3: 13
+DEAL:0:all::n_active_cells on level=4: 15
+DEAL:0:all::n_active_cells on level=5: 23
+DEAL:0:all::n_active_cells on level=6: 372
+DEAL:0:all::
+DEAL:0:all::n_dofs: 4813
+DEAL:0:all::n_locally_owned_dofs: 1233
+DEAL:0:all::n_dofs on level=0: 9
+DEAL:0:all::n_dofs on level=1: 25
+DEAL:0:all::n_dofs on level=2: 81
+DEAL:0:all::n_dofs on level=3: 289
+DEAL:0:all::n_dofs on level=4: 441
+DEAL:0:all::n_dofs on level=5: 1369
+DEAL:0:all::n_dofs on level=6: 4225
+DEAL:0:all::n_locally_owned_mg_dofs on level=0: 9
+DEAL:0:all::n_locally_owned_mg_dofs on level=1: 9
+DEAL:0:all::n_locally_owned_mg_dofs on level=2: 15
+DEAL:0:all::n_locally_owned_mg_dofs on level=3: 31
+DEAL:0:all::n_locally_owned_mg_dofs on level=4: 93
+DEAL:0:all::n_locally_owned_mg_dofs on level=5: 329
+DEAL:0:all::n_locally_owned_mg_dofs on level=6: 1233
+DEAL:0:all::
+DEAL:0:all::n_dofs: 4813
+DEAL:0:all::n_locally_owned_dofs: 1233
+DEAL:0:all::n_dofs on level=0: 9
+DEAL:0:all::n_dofs on level=1: 25
+DEAL:0:all::n_dofs on level=2: 81
+DEAL:0:all::n_dofs on level=3: 289
+DEAL:0:all::n_dofs on level=4: 441
+DEAL:0:all::n_dofs on level=5: 1369
+DEAL:0:all::n_dofs on level=6: 4225
+DEAL:0:all::n_locally_owned_mg_dofs on level=0: 9
+DEAL:0:all::n_locally_owned_mg_dofs on level=1: 9
+DEAL:0:all::n_locally_owned_mg_dofs on level=2: 15
+DEAL:0:all::n_locally_owned_mg_dofs on level=3: 31
+DEAL:0:all::n_locally_owned_mg_dofs on level=4: 93
+DEAL:0:all::n_locally_owned_mg_dofs on level=5: 329
+DEAL:0:all::n_locally_owned_mg_dofs on level=6: 1233
+DEAL:0:all::
+
+DEAL:1:all::n_levels: 7
+DEAL:1:all::n_cells: 713
+DEAL:1:all::n_active_cells: 535
+DEAL:1:all::n_cells on level=0: 1
+DEAL:1:all::n_cells on level=1: 4
+DEAL:1:all::n_cells on level=2: 16
+DEAL:1:all::n_cells on level=3: 32
+DEAL:1:all::n_cells on level=4: 64
+DEAL:1:all::n_cells on level=5: 168
+DEAL:1:all::n_cells on level=6: 428
+DEAL:1:all::n_active_cells on level=0: 0
+DEAL:1:all::n_active_cells on level=1: 0
+DEAL:1:all::n_active_cells on level=2: 8
+DEAL:1:all::n_active_cells on level=3: 16
+DEAL:1:all::n_active_cells on level=4: 22
+DEAL:1:all::n_active_cells on level=5: 61
+DEAL:1:all::n_active_cells on level=6: 428
+DEAL:1:all::
+DEAL:1:all::n_levels: 7
+DEAL:1:all::n_cells: 713
+DEAL:1:all::n_active_cells: 535
+DEAL:1:all::n_cells on level=0: 1
+DEAL:1:all::n_cells on level=1: 4
+DEAL:1:all::n_cells on level=2: 16
+DEAL:1:all::n_cells on level=3: 32
+DEAL:1:all::n_cells on level=4: 64
+DEAL:1:all::n_cells on level=5: 168
+DEAL:1:all::n_cells on level=6: 428
+DEAL:1:all::n_active_cells on level=0: 0
+DEAL:1:all::n_active_cells on level=1: 0
+DEAL:1:all::n_active_cells on level=2: 8
+DEAL:1:all::n_active_cells on level=3: 16
+DEAL:1:all::n_active_cells on level=4: 22
+DEAL:1:all::n_active_cells on level=5: 61
+DEAL:1:all::n_active_cells on level=6: 428
+DEAL:1:all::
+DEAL:1:all::n_dofs: 4813
+DEAL:1:all::n_locally_owned_dofs: 1184
+DEAL:1:all::n_dofs on level=0: 9
+DEAL:1:all::n_dofs on level=1: 25
+DEAL:1:all::n_dofs on level=2: 81
+DEAL:1:all::n_dofs on level=3: 289
+DEAL:1:all::n_dofs on level=4: 441
+DEAL:1:all::n_dofs on level=5: 1369
+DEAL:1:all::n_dofs on level=6: 4225
+DEAL:1:all::n_locally_owned_mg_dofs on level=0: 0
+DEAL:1:all::n_locally_owned_mg_dofs on level=1: 0
+DEAL:1:all::n_locally_owned_mg_dofs on level=2: 6
+DEAL:1:all::n_locally_owned_mg_dofs on level=3: 20
+DEAL:1:all::n_locally_owned_mg_dofs on level=4: 80
+DEAL:1:all::n_locally_owned_mg_dofs on level=5: 304
+DEAL:1:all::n_locally_owned_mg_dofs on level=6: 1184
+DEAL:1:all::
+DEAL:1:all::n_dofs: 4813
+DEAL:1:all::n_locally_owned_dofs: 1184
+DEAL:1:all::n_dofs on level=0: 9
+DEAL:1:all::n_dofs on level=1: 25
+DEAL:1:all::n_dofs on level=2: 81
+DEAL:1:all::n_dofs on level=3: 289
+DEAL:1:all::n_dofs on level=4: 441
+DEAL:1:all::n_dofs on level=5: 1369
+DEAL:1:all::n_dofs on level=6: 4225
+DEAL:1:all::n_locally_owned_mg_dofs on level=0: 0
+DEAL:1:all::n_locally_owned_mg_dofs on level=1: 0
+DEAL:1:all::n_locally_owned_mg_dofs on level=2: 6
+DEAL:1:all::n_locally_owned_mg_dofs on level=3: 20
+DEAL:1:all::n_locally_owned_mg_dofs on level=4: 80
+DEAL:1:all::n_locally_owned_mg_dofs on level=5: 304
+DEAL:1:all::n_locally_owned_mg_dofs on level=6: 1184
+DEAL:1:all::
+
+
+DEAL:2:all::n_levels: 7
+DEAL:2:all::n_cells: 717
+DEAL:2:all::n_active_cells: 538
+DEAL:2:all::n_cells on level=0: 1
+DEAL:2:all::n_cells on level=1: 4
+DEAL:2:all::n_cells on level=2: 16
+DEAL:2:all::n_cells on level=3: 36
+DEAL:2:all::n_cells on level=4: 72
+DEAL:2:all::n_cells on level=5: 168
+DEAL:2:all::n_cells on level=6: 420
+DEAL:2:all::n_active_cells on level=0: 0
+DEAL:2:all::n_active_cells on level=1: 0
+DEAL:2:all::n_active_cells on level=2: 7
+DEAL:2:all::n_active_cells on level=3: 18
+DEAL:2:all::n_active_cells on level=4: 30
+DEAL:2:all::n_active_cells on level=5: 63
+DEAL:2:all::n_active_cells on level=6: 420
+DEAL:2:all::
+DEAL:2:all::n_levels: 7
+DEAL:2:all::n_cells: 717
+DEAL:2:all::n_active_cells: 538
+DEAL:2:all::n_cells on level=0: 1
+DEAL:2:all::n_cells on level=1: 4
+DEAL:2:all::n_cells on level=2: 16
+DEAL:2:all::n_cells on level=3: 36
+DEAL:2:all::n_cells on level=4: 72
+DEAL:2:all::n_cells on level=5: 168
+DEAL:2:all::n_cells on level=6: 420
+DEAL:2:all::n_active_cells on level=0: 0
+DEAL:2:all::n_active_cells on level=1: 0
+DEAL:2:all::n_active_cells on level=2: 7
+DEAL:2:all::n_active_cells on level=3: 18
+DEAL:2:all::n_active_cells on level=4: 30
+DEAL:2:all::n_active_cells on level=5: 63
+DEAL:2:all::n_active_cells on level=6: 420
+DEAL:2:all::
+DEAL:2:all::n_dofs: 4813
+DEAL:2:all::n_locally_owned_dofs: 1168
+DEAL:2:all::n_dofs on level=0: 9
+DEAL:2:all::n_dofs on level=1: 25
+DEAL:2:all::n_dofs on level=2: 81
+DEAL:2:all::n_dofs on level=3: 289
+DEAL:2:all::n_dofs on level=4: 441
+DEAL:2:all::n_dofs on level=5: 1369
+DEAL:2:all::n_dofs on level=6: 4225
+DEAL:2:all::n_locally_owned_mg_dofs on level=0: 0
+DEAL:2:all::n_locally_owned_mg_dofs on level=1: 0
+DEAL:2:all::n_locally_owned_mg_dofs on level=2: 4
+DEAL:2:all::n_locally_owned_mg_dofs on level=3: 22
+DEAL:2:all::n_locally_owned_mg_dofs on level=4: 76
+DEAL:2:all::n_locally_owned_mg_dofs on level=5: 296
+DEAL:2:all::n_locally_owned_mg_dofs on level=6: 1168
+DEAL:2:all::
+DEAL:2:all::n_dofs: 4813
+DEAL:2:all::n_locally_owned_dofs: 1168
+DEAL:2:all::n_dofs on level=0: 9
+DEAL:2:all::n_dofs on level=1: 25
+DEAL:2:all::n_dofs on level=2: 81
+DEAL:2:all::n_dofs on level=3: 289
+DEAL:2:all::n_dofs on level=4: 441
+DEAL:2:all::n_dofs on level=5: 1369
+DEAL:2:all::n_dofs on level=6: 4225
+DEAL:2:all::n_locally_owned_mg_dofs on level=0: 0
+DEAL:2:all::n_locally_owned_mg_dofs on level=1: 0
+DEAL:2:all::n_locally_owned_mg_dofs on level=2: 4
+DEAL:2:all::n_locally_owned_mg_dofs on level=3: 22
+DEAL:2:all::n_locally_owned_mg_dofs on level=4: 76
+DEAL:2:all::n_locally_owned_mg_dofs on level=5: 296
+DEAL:2:all::n_locally_owned_mg_dofs on level=6: 1168
+DEAL:2:all::
+
+
+DEAL:3:all::n_levels: 7
+DEAL:3:all::n_cells: 613
+DEAL:3:all::n_active_cells: 460
+DEAL:3:all::n_cells on level=0: 1
+DEAL:3:all::n_cells on level=1: 4
+DEAL:3:all::n_cells on level=2: 16
+DEAL:3:all::n_cells on level=3: 64
+DEAL:3:all::n_cells on level=4: 76
+DEAL:3:all::n_cells on level=5: 168
+DEAL:3:all::n_cells on level=6: 284
+DEAL:3:all::n_active_cells on level=0: 0
+DEAL:3:all::n_active_cells on level=1: 0
+DEAL:3:all::n_active_cells on level=2: 0
+DEAL:3:all::n_active_cells on level=3: 45
+DEAL:3:all::n_active_cells on level=4: 34
+DEAL:3:all::n_active_cells on level=5: 97
+DEAL:3:all::n_active_cells on level=6: 284
+DEAL:3:all::
+DEAL:3:all::n_levels: 7
+DEAL:3:all::n_cells: 613
+DEAL:3:all::n_active_cells: 460
+DEAL:3:all::n_cells on level=0: 1
+DEAL:3:all::n_cells on level=1: 4
+DEAL:3:all::n_cells on level=2: 16
+DEAL:3:all::n_cells on level=3: 64
+DEAL:3:all::n_cells on level=4: 76
+DEAL:3:all::n_cells on level=5: 168
+DEAL:3:all::n_cells on level=6: 284
+DEAL:3:all::n_active_cells on level=0: 0
+DEAL:3:all::n_active_cells on level=1: 0
+DEAL:3:all::n_active_cells on level=2: 0
+DEAL:3:all::n_active_cells on level=3: 45
+DEAL:3:all::n_active_cells on level=4: 34
+DEAL:3:all::n_active_cells on level=5: 97
+DEAL:3:all::n_active_cells on level=6: 284
+DEAL:3:all::
+DEAL:3:all::n_dofs: 4813
+DEAL:3:all::n_locally_owned_dofs: 1228
+DEAL:3:all::n_dofs on level=0: 9
+DEAL:3:all::n_dofs on level=1: 25
+DEAL:3:all::n_dofs on level=2: 81
+DEAL:3:all::n_dofs on level=3: 289
+DEAL:3:all::n_dofs on level=4: 441
+DEAL:3:all::n_dofs on level=5: 1369
+DEAL:3:all::n_dofs on level=6: 4225
+DEAL:3:all::n_locally_owned_mg_dofs on level=0: 0
+DEAL:3:all::n_locally_owned_mg_dofs on level=1: 16
+DEAL:3:all::n_locally_owned_mg_dofs on level=2: 56
+DEAL:3:all::n_locally_owned_mg_dofs on level=3: 216
+DEAL:3:all::n_locally_owned_mg_dofs on level=4: 192
+DEAL:3:all::n_locally_owned_mg_dofs on level=5: 440
+DEAL:3:all::n_locally_owned_mg_dofs on level=6: 640
+DEAL:3:all::
+DEAL:3:all::n_dofs: 4813
+DEAL:3:all::n_locally_owned_dofs: 1228
+DEAL:3:all::n_dofs on level=0: 9
+DEAL:3:all::n_dofs on level=1: 25
+DEAL:3:all::n_dofs on level=2: 81
+DEAL:3:all::n_dofs on level=3: 289
+DEAL:3:all::n_dofs on level=4: 441
+DEAL:3:all::n_dofs on level=5: 1369
+DEAL:3:all::n_dofs on level=6: 4225
+DEAL:3:all::n_locally_owned_mg_dofs on level=0: 0
+DEAL:3:all::n_locally_owned_mg_dofs on level=1: 16
+DEAL:3:all::n_locally_owned_mg_dofs on level=2: 56
+DEAL:3:all::n_locally_owned_mg_dofs on level=3: 216
+DEAL:3:all::n_locally_owned_mg_dofs on level=4: 192
+DEAL:3:all::n_locally_owned_mg_dofs on level=5: 440
+DEAL:3:all::n_locally_owned_mg_dofs on level=6: 640
+DEAL:3:all::
+