]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Fix sorting of cells in p:f:t 10092/head
authorPeter Munch <peterrmuench@gmail.com>
Fri, 8 May 2020 15:33:30 +0000 (17:33 +0200)
committerPeter Munch <peterrmuench@gmail.com>
Fri, 8 May 2020 16:14:06 +0000 (18:14 +0200)
source/distributed/fully_distributed_tria.cc
tests/fullydistributed_grids/copy_serial_tria_10.cc [new file with mode: 0644]
tests/fullydistributed_grids/copy_serial_tria_10.mpirun=1.output [new file with mode: 0644]

index 5f5e90d42ae2fbae9c6d96f3f2e377092a314b0a..23a393d1f854992c4460ffd79816d3fd53a04869 100644 (file)
@@ -140,6 +140,34 @@ namespace parallel
           // create a copy of cell_infos such that we can sort them
           auto cell_infos = construction_data.cell_infos;
 
+          // sort cell_infos on each level separately (as done in
+          // dealii::Triangulation::create_triangulation())
+          for (auto &cell_info : cell_infos)
+            std::sort(cell_info.begin(),
+                      cell_info.end(),
+                      [&](TriangulationDescription::CellData<dim> a,
+                          TriangulationDescription::CellData<dim> b) {
+                        const CellId a_id(a.id);
+                        const CellId b_id(b.id);
+
+                        const auto a_coarse_cell_index =
+                          this->coarse_cell_id_to_coarse_cell_index(
+                            a_id.get_coarse_cell_id());
+                        const auto b_coarse_cell_index =
+                          this->coarse_cell_id_to_coarse_cell_index(
+                            b_id.get_coarse_cell_id());
+
+                        // according to their coarse-cell index and if that is
+                        // same according to their cell id (the result is that
+                        // cells on each level are sorted according to their
+                        // index on that level - what we need in the following
+                        // operations)
+                        if (a_coarse_cell_index != b_coarse_cell_index)
+                          return a_coarse_cell_index < b_coarse_cell_index;
+                        else
+                          return a_id < b_id;
+                      });
+
           // 4a) set all cells artificial (and set the actual
           //     (level_)subdomain_ids in the next step)
           for (auto cell = this->begin(); cell != this->end(); ++cell)
diff --git a/tests/fullydistributed_grids/copy_serial_tria_10.cc b/tests/fullydistributed_grids/copy_serial_tria_10.cc
new file mode 100644 (file)
index 0000000..91114c7
--- /dev/null
@@ -0,0 +1,102 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2019 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.
+//
+// ---------------------------------------------------------------------
+
+
+// Similar as copy_serial_tria_01 but permute the cells on the coarse grid
+// s.t. the cells in the auto-created TriangulationDescription::cell_infos
+// does not match the order of coarse cells and sorting is required in
+// Triangulation::create_triangulation() (see also PR #10092).
+
+#include <deal.II/base/mpi.h>
+
+#include <deal.II/distributed/fully_distributed_tria.h>
+#include <deal.II/distributed/shared_tria.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.h>
+#include <deal.II/grid/tria_description.h>
+
+#include "./tests.h"
+
+using namespace dealii;
+
+template <int dim>
+void
+test(int n_refinements, MPI_Comm comm)
+{
+  // create serial triangulation
+  Triangulation<dim> basetria;
+  GridGenerator::hyper_L(basetria);
+  basetria.refine_global(n_refinements);
+
+  GridTools::partition_triangulation_zorder(
+    Utilities::MPI::n_mpi_processes(comm), basetria);
+
+  // create instance of pft
+  parallel::fullydistributed::Triangulation<dim> tria_pft(comm);
+
+  // extract relevant information form serial triangulation
+  auto construction_data =
+    TriangulationDescription::Utilities::create_description_from_triangulation(
+      basetria, comm);
+
+  for (unsigned int i = 0; i < construction_data.coarse_cells.size() / 2; ++i)
+    {
+      std::swap(construction_data.coarse_cells[i],
+                construction_data
+                  .coarse_cells[construction_data.coarse_cells.size() - 1 - i]);
+      std::swap(construction_data.coarse_cell_index_to_coarse_cell_id[i],
+                construction_data.coarse_cell_index_to_coarse_cell_id
+                  [construction_data.coarse_cells.size() - 1 - i]);
+    }
+
+  // actually create triangulation
+  tria_pft.create_triangulation(construction_data);
+
+  // test triangulation
+  FE_Q<dim>       fe(2);
+  DoFHandler<dim> dof_handler(tria_pft);
+  dof_handler.distribute_dofs(fe);
+
+  // print statistics
+  print_statistics(tria_pft);
+  print_statistics(dof_handler);
+}
+
+int
+main(int argc, char *argv[])
+{
+  Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);
+  MPILogInitAll                    all;
+
+  const int      n_refinements = 3;
+  const MPI_Comm comm          = MPI_COMM_WORLD;
+
+  {
+    deallog.push("2d");
+    test<2>(n_refinements, comm);
+    deallog.pop();
+  }
+  {
+    deallog.push("3d");
+    test<3>(n_refinements, comm);
+    deallog.pop();
+  }
+}
diff --git a/tests/fullydistributed_grids/copy_serial_tria_10.mpirun=1.output b/tests/fullydistributed_grids/copy_serial_tria_10.mpirun=1.output
new file mode 100644 (file)
index 0000000..cb1dcf4
--- /dev/null
@@ -0,0 +1,15 @@
+
+DEAL:0:2d::n_levels:                  4
+DEAL:0:2d::n_cells:                   255
+DEAL:0:2d::n_active_cells:            192
+DEAL:0:2d::
+DEAL:0:2d::n_dofs:                             833
+DEAL:0:2d::n_locally_owned_dofs:               833
+DEAL:0:2d::
+DEAL:0:3d::n_levels:                  4
+DEAL:0:3d::n_cells:                   4095
+DEAL:0:3d::n_active_cells:            3584
+DEAL:0:3d::
+DEAL:0:3d::n_dofs:                             31841
+DEAL:0:3d::n_locally_owned_dofs:               31841
+DEAL:0:3d::

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.