]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Coninue
authorPeter Munch <peterrmuench@gmail.com>
Tue, 7 Dec 2021 20:58:42 +0000 (21:58 +0100)
committerPeter Munch <peterrmuench@gmail.com>
Tue, 7 Dec 2021 21:03:54 +0000 (22:03 +0100)
source/distributed/fully_distributed_tria.cc
source/grid/tria.cc
source/grid/tria_description.cc
tests/fullydistributed_grids/repartitioning_08.cc [new file with mode: 0644]
tests/fullydistributed_grids/repartitioning_08.with_p4est=true.mpirun=4.output [new file with mode: 0644]

index 2a7496d5ec07eac3744ee6365c46867174f6c758..7e6df2b8158184e5ac01601fb6c294be171734e5 100644 (file)
@@ -335,6 +335,10 @@ namespace parallel
           this->partitioner_distributed->partition(*this),
           this->settings);
 
+      this->clear();
+      this->coarse_cell_id_to_coarse_cell_index_vector.clear();
+      this->coarse_cell_index_to_coarse_cell_id_vector.clear();
+
       this->create_triangulation(construction_data);
 
       this->signals.post_distributed_repartition();
index 57fce4ad707c3fcc9005d3036b51d6ba8cc1d849..76b33db228cf4480faa387a3d23d6c7dd66dee44 100644 (file)
@@ -11871,11 +11871,8 @@ Triangulation<dim, spacedim>::create_triangulation(
 
           // boundary ids
           for (auto pair : cell_info->boundary_ids)
-            {
-              Assert(cell->at_boundary(pair.first),
-                     ExcMessage("Cell face is not on the boundary!"));
+            if (cell->face(pair.first)->at_boundary())
               cell->face(pair.first)->set_boundary_id(pair.second);
-            }
         }
     }
 }
index 1606c0df6481e7451510f013e7c93acb132d8aa3..e8c5041c0b514d2eccc3fb54207f35bb3acd7b1a 100644 (file)
@@ -163,9 +163,6 @@ namespace TriangulationDescription
               const auto comp = [](const auto &a, const auto &b) {
                 const double tolerance = 1e-10;
 
-                if (a.distance(b) < tolerance)
-                  return false;
-
                 for (unsigned int i = 0; i < spacedim; ++i)
                   if (std::abs(a[i] - b[i]) > tolerance)
                     return a[i] < b[i];
@@ -191,15 +188,17 @@ namespace TriangulationDescription
                 else
                   map_i[p.first] = map[p.second];
 
-              for (auto cell : other.coarse_cells)
+              auto other_coarse_cells_copy = other.coarse_cells;
+
+              for (auto &cell : other_coarse_cells_copy)
                 {
                   for (auto &v : cell.vertices)
                     v = map_i[v];
                 }
 
               this->coarse_cells.insert(this->coarse_cells.end(),
-                                        other.coarse_cells.begin(),
-                                        other.coarse_cells.end());
+                                        other_coarse_cells_copy.begin(),
+                                        other_coarse_cells_copy.end());
             }
           else
             {
@@ -280,7 +279,13 @@ namespace TriangulationDescription
               std::unique(this->coarse_cell_vertices.begin(),
                           this->coarse_cell_vertices.end(),
                           [](const auto &a, const auto &b) {
-                            return a.first == b.first;
+                            if (a.first == b.first)
+                              {
+                                Assert(a.second.distance(b.second) < 10e-8,
+                                       ExcInternalError());
+                                return true;
+                              }
+                            return false;
                           }),
               this->coarse_cell_vertices.end());
           }
@@ -1083,7 +1088,7 @@ namespace TriangulationDescription
         partition.get_mpi_communicator(),
         dynamic_cast<
           const parallel::fullydistributed::Triangulation<dim, spacedim> *>(
-          &tria) != nullptr);
+          &tria) == nullptr);
 
       // remove redundant entries
       description_merged.reduce();
diff --git a/tests/fullydistributed_grids/repartitioning_08.cc b/tests/fullydistributed_grids/repartitioning_08.cc
new file mode 100644 (file)
index 0000000..5905486
--- /dev/null
@@ -0,0 +1,153 @@
+// ---------------------------------------------------------------------
+//
+// 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()
+// with repartitioning capabilities.
+
+#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;
+
+
+template <int dim, int spacedim = dim>
+class MyPolicy : public RepartitioningPolicyTools::Base<dim, spacedim>
+{
+public:
+  MyPolicy(const MPI_Comm &comm, const unsigned int direction)
+    : comm(comm)
+    , direction(direction)
+  {}
+
+  virtual LinearAlgebra::distributed::Vector<double>
+  partition(const Triangulation<dim, spacedim> &tria_in) const override
+  {
+    const auto tria =
+      dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(
+        &tria_in);
+
+    Assert(tria, ExcNotImplemented());
+
+    LinearAlgebra::distributed::Vector<double> partition(
+      tria->global_active_cell_index_partitioner().lock());
+
+    const unsigned int n_partitions = Utilities::MPI::n_mpi_processes(comm);
+
+    for (const auto &cell :
+         tria_in.active_cell_iterators() | IteratorFilters::LocallyOwnedCell())
+      partition[cell->global_active_cell_index()] =
+        std::floor(cell->center()[direction] * n_partitions);
+
+    partition.update_ghost_values();
+
+    return partition;
+  }
+
+private:
+  const MPI_Comm &   comm;
+  const unsigned int direction;
+};
+
+
+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);
+  GridGenerator::subdivided_hyper_cube(tria, 4);
+  tria.refine_global(2);
+
+  MyPolicy<dim> policy_0(comm, 0);
+  MyPolicy<dim> policy_1(comm, 1);
+
+  const auto partition_0 = policy_0.partition(tria);
+
+  const auto settings =
+    TriangulationDescription::Settings::construct_multigrid_hierarchy;
+
+  // repartition triangulation so that it has strided partitioning
+  const auto construction_data =
+    TriangulationDescription::Utilities::create_description_from_triangulation(
+      tria, partition_0, settings);
+
+  parallel::fullydistributed::Triangulation<dim> tria_pft(comm);
+  tria_pft.create_triangulation(construction_data);
+
+  {
+    FE_Q<dim>       fe(2);
+    DoFHandler<dim> dof_handler(tria_pft);
+    dof_handler.distribute_dofs(fe);
+    dof_handler.distribute_mg_dofs();
+
+    // print statistics
+    print_statistics(tria_pft);
+    print_statistics(dof_handler);
+  }
+
+  GridOut go;
+  go.write_mesh_per_processor_as_vtu(tria_pft, "mesh_old");
+
+  tria_pft.set_partitioner(policy_1, settings);
+
+  tria_pft.repartition();
+
+  go.write_mesh_per_processor_as_vtu(tria_pft, "mesh_new");
+
+  {
+    FE_Q<dim>       fe(2);
+    DoFHandler<dim> dof_handler(tria_pft);
+    dof_handler.distribute_dofs(fe);
+    dof_handler.distribute_mg_dofs();
+
+    // print statistics
+    print_statistics(tria_pft);
+    print_statistics(dof_handler);
+  }
+}
+
+
+
+int
+main(int argc, char **argv)
+{
+  Utilities::MPI::MPI_InitFinalize mpi(argc, argv, 1);
+  MPILogInitAll                    all;
+
+  MPI_Comm comm = MPI_COMM_WORLD;
+
+  test<2>(comm);
+  deallog.pop();
+}
diff --git a/tests/fullydistributed_grids/repartitioning_08.with_p4est=true.mpirun=4.output b/tests/fullydistributed_grids/repartitioning_08.with_p4est=true.mpirun=4.output
new file mode 100644 (file)
index 0000000..658d881
--- /dev/null
@@ -0,0 +1,63 @@
+
+DEAL:0::n_levels:                  3
+DEAL:0::n_cells:                   136
+DEAL:0::n_active_cells:            104
+DEAL:0::
+DEAL:0::n_dofs:                             1089
+DEAL:0::n_locally_owned_dofs:               297
+DEAL:0::
+DEAL:0::n_levels:                  3
+DEAL:0::n_cells:                   136
+DEAL:0::n_active_cells:            104
+DEAL:0::
+DEAL:0::n_dofs:                             1089
+DEAL:0::n_locally_owned_dofs:               297
+DEAL:0::
+
+DEAL:1::n_levels:                  3
+DEAL:1::n_cells:                   188
+DEAL:1::n_active_cells:            144
+DEAL:1::
+DEAL:1::n_dofs:                             1089
+DEAL:1::n_locally_owned_dofs:               264
+DEAL:1::
+DEAL:1::n_levels:                  3
+DEAL:1::n_cells:                   188
+DEAL:1::n_active_cells:            144
+DEAL:1::
+DEAL:1::n_dofs:                             1089
+DEAL:1::n_locally_owned_dofs:               264
+DEAL:1::
+
+
+DEAL:2::n_levels:                  3
+DEAL:2::n_cells:                   188
+DEAL:2::n_active_cells:            144
+DEAL:2::
+DEAL:2::n_dofs:                             1089
+DEAL:2::n_locally_owned_dofs:               264
+DEAL:2::
+DEAL:2::n_levels:                  3
+DEAL:2::n_cells:                   188
+DEAL:2::n_active_cells:            144
+DEAL:2::
+DEAL:2::n_dofs:                             1089
+DEAL:2::n_locally_owned_dofs:               264
+DEAL:2::
+
+
+DEAL:3::n_levels:                  3
+DEAL:3::n_cells:                   136
+DEAL:3::n_active_cells:            104
+DEAL:3::
+DEAL:3::n_dofs:                             1089
+DEAL:3::n_locally_owned_dofs:               264
+DEAL:3::
+DEAL:3::n_levels:                  3
+DEAL:3::n_cells:                   136
+DEAL:3::n_active_cells:            104
+DEAL:3::
+DEAL:3::n_dofs:                             1089
+DEAL:3::n_locally_owned_dofs:               264
+DEAL:3::
+

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.