]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Fix treatment of periodic faces in communicate_locally_moved_vertices 5612/head
authorDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Mon, 11 Dec 2017 14:53:52 +0000 (15:53 +0100)
committerDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Wed, 13 Dec 2017 08:44:35 +0000 (09:44 +0100)
doc/news/changes/minor/20171211DanielArndtSambitDas [new file with mode: 0644]
source/distributed/tria.cc
tests/mpi/communicate_moved_vertices_03.cc [new file with mode: 0644]
tests/mpi/communicate_moved_vertices_03.mpirun=1.output [new file with mode: 0644]
tests/mpi/communicate_moved_vertices_03.mpirun=3.output [new file with mode: 0644]

diff --git a/doc/news/changes/minor/20171211DanielArndtSambitDas b/doc/news/changes/minor/20171211DanielArndtSambitDas
new file mode 100644 (file)
index 0000000..4926015
--- /dev/null
@@ -0,0 +1,4 @@
+Fixed: parallel::distributed::Triangulation::communicate_locally_moved_vertices
+treats periodic faces correctly now.
+<br>
+(Daniel Arndt, Sambit Das, 2017/12/11)
index eaffda540ff97f152cc325681a7cdd4e275fd364..2f6731e24d2c9b1a0720fb22a36fbedab6c44454 100644 (file)
@@ -3064,23 +3064,14 @@ namespace parallel
       }
 #endif
 
-      std::map<unsigned int, std::set<dealii::types::subdomain_id> >
-      vertices_with_ghost_neighbors;
-
       // First find out which process should receive which vertices.
-      // these are specifically the ones that sit on ghost cells and,
-      // among these, the ones that we own locally
-      for (typename Triangulation<dim,spacedim>::active_cell_iterator
-           cell=this->begin_active(); cell!=this->end();
-           ++cell)
-        if (cell->is_ghost())
-          for (unsigned int vertex_no=0;
-               vertex_no<GeometryInfo<dim>::vertices_per_cell; ++vertex_no)
-            {
-              const unsigned int process_local_vertex_no = cell->vertex_index(vertex_no);
-              vertices_with_ghost_neighbors[process_local_vertex_no].insert
-              (cell->subdomain_id());
-            }
+      // These are specifically the ones that are located on cells at the
+      // boundary of the subdomain this process owns and the receiving
+      // process taking periodic faces into account.
+      // Here, it is sufficient to collect all vertices that are located
+      // at that boundary.
+      const std::map<unsigned int, std::set<dealii::types::subdomain_id> >
+      vertices_with_ghost_neighbors = compute_vertices_with_ghost_neighbors ();
 
       // now collect cells and their vertices
       // for the interested neighbors
diff --git a/tests/mpi/communicate_moved_vertices_03.cc b/tests/mpi/communicate_moved_vertices_03.cc
new file mode 100644 (file)
index 0000000..5a8bfda
--- /dev/null
@@ -0,0 +1,99 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2017 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 at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+// test Triangulation::communicate_locally_moved_vertices() on a
+// Triangulation with periodic faces.
+
+#include "../tests.h"
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/grid_tools.h>
+#include <deal.II/grid/grid_out.h>
+#include <deal.II/distributed/tria.h>
+
+template <int dim>
+void test()
+{
+  deallog << "Testing " << Utilities::int_to_string(dim,1) << "D" << std::endl;
+
+  parallel::distributed::Triangulation<dim> triangulation(MPI_COMM_WORLD);
+  GridGenerator::hyper_cube (triangulation, -1., 1., true);
+
+  std::vector<GridTools::PeriodicFacePair<typename parallel::distributed::Triangulation<dim>::cell_iterator> > periodicity_vector;
+  for (int i = 0; i < dim; ++i)
+    GridTools::collect_periodic_faces(triangulation, 2*i, 2*i+1, i, periodicity_vector);
+  triangulation.add_periodicity(periodicity_vector);
+
+  triangulation.refine_global(3);
+
+  auto cell = triangulation.begin_active();
+  const auto endc = triangulation.end();
+
+  std::map<unsigned int, Point<dim>> non_artificial_vertices_old;
+  for (cell = triangulation.begin_active(); cell!=endc; ++cell)
+    if (!cell->is_artificial())
+      for (unsigned int vertex_no=0; vertex_no<GeometryInfo<dim>::vertices_per_cell; ++vertex_no)
+        non_artificial_vertices_old[cell->vertex_index(vertex_no)] = cell->vertex(vertex_no);
+
+  std::vector<bool> vertex_moved(triangulation.n_vertices(), false);
+  const std::vector<bool> locally_owned_vertices = GridTools::get_locally_owned_vertices(triangulation);
+  for (cell = triangulation.begin_active(); cell!=endc; ++cell)
+    if (cell->is_locally_owned())
+      for (unsigned int vertex_no=0; vertex_no<GeometryInfo<dim>::vertices_per_cell;
+           ++vertex_no)
+        {
+          const unsigned global_vertex_no = cell->vertex_index(vertex_no);
+          if (!vertex_moved[global_vertex_no] && locally_owned_vertices[global_vertex_no])
+            {
+              cell->vertex(vertex_no)(0) += 1.e-1;
+              vertex_moved[global_vertex_no] = true;
+            }
+        }
+
+  triangulation.communicate_locally_moved_vertices(vertex_moved);
+
+  std::map<unsigned int, Point<dim>> non_artificial_vertices_new;
+  for (cell = triangulation.begin_active(); cell!=endc; ++cell)
+    if (!cell->is_artificial())
+      for (unsigned int vertex_no=0; vertex_no<GeometryInfo<dim>::vertices_per_cell;
+           ++vertex_no)
+        {
+          Point<dim> point = cell->vertex(vertex_no);
+          point(0)-=1.e-1;
+          non_artificial_vertices_new[cell->vertex_index(vertex_no)] = point;
+        }
+
+  for (const auto &pair: non_artificial_vertices_new)
+    if ((non_artificial_vertices_old[pair.first]-pair.second).norm_square()>1.e-6)
+      {
+        deallog << pair.first << ": "
+                << non_artificial_vertices_old[pair.first] << " vs. " << pair.second << std::endl;
+        AssertThrow(false, ExcMessage("Some of the vertices on ghost cell were not moved correctly!"));
+      }
+
+  /*std::string filename("grid"+Utilities::int_to_string(Utilities::MPI::this_mpi_process(MPI_COMM_WORLD),2)+".vtu");
+  std::ofstream out(filename.c_str());
+  GridOut grid_out;
+  grid_out.write_vtu(triangulation, out);*/
+  deallog << Utilities::int_to_string(dim,1) << "D OK" << std::endl;
+}
+
+int main (int argc, char *argv[])
+{
+  Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv);
+  MPILogInitAll log;
+  test<2>();
+  test<3>();
+  return 0;
+}
diff --git a/tests/mpi/communicate_moved_vertices_03.mpirun=1.output b/tests/mpi/communicate_moved_vertices_03.mpirun=1.output
new file mode 100644 (file)
index 0000000..47c8c9e
--- /dev/null
@@ -0,0 +1,5 @@
+
+DEAL:0::Testing 2D
+DEAL:0::2D OK
+DEAL:0::Testing 3D
+DEAL:0::3D OK
diff --git a/tests/mpi/communicate_moved_vertices_03.mpirun=3.output b/tests/mpi/communicate_moved_vertices_03.mpirun=3.output
new file mode 100644 (file)
index 0000000..020bc46
--- /dev/null
@@ -0,0 +1,17 @@
+
+DEAL:0::Testing 2D
+DEAL:0::2D OK
+DEAL:0::Testing 3D
+DEAL:0::3D OK
+
+DEAL:1::Testing 2D
+DEAL:1::2D OK
+DEAL:1::Testing 3D
+DEAL:1::3D OK
+
+
+DEAL:2::Testing 2D
+DEAL:2::2D OK
+DEAL:2::Testing 3D
+DEAL:2::3D OK
+

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.