]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Fix a segmentation fault in DoFHandler::locally_owned_dofs(). 2388/head
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Sun, 20 Mar 2016 18:52:31 +0000 (13:52 -0500)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Sun, 20 Mar 2016 18:52:31 +0000 (13:52 -0500)
This was caused by an incorrect size computation inside DoFTools::locally_owned_dofs_per_subdomain().

doc/news/changes.h
source/dofs/dof_tools.cc
tests/mpi/locally_owned_dofs_per_subdomain.cc [new file with mode: 0644]
tests/mpi/locally_owned_dofs_per_subdomain.mpirun=2.output [new file with mode: 0644]

index ec5a260aee1a2818124f0205d8f8e25c497899d3..f4384b45a23c80d2115cc9141a676db2d795638d 100644 (file)
@@ -104,6 +104,13 @@ inconvenience this causes.
 <h3>Specific improvements</h3>
 
 <ol>
+ <li> Fixed: DoFHandler::locally_owned_dofs() could create a segmentation
+ fault in cases where some processors do not own any cells. This was caused
+ by an incorrect computation in DoFTools::locally_owned_dofs_per_subdomain().
+ <br>
+ (Wolfgang Bangerth, 2016/03/20)
+ </li>
+
  <li> Improved: The distribution of degrees of freedom on multigrid levels,
  DoFHandler::distribute_mg_dofs(), contained a few steps that scaled
  quadratically in the number of local cells for certain configurations. These
index 1aa3e161a3b944da1f5a3f64f4f976356354f2ba..502e791cb585c2fd95e800437eed97447bd92752 100644 (file)
@@ -1145,10 +1145,22 @@ namespace DoFTools
     std::vector< dealii::types::subdomain_id > subdomain_association (dof_handler.n_dofs ());
     dealii::DoFTools::get_subdomain_association (dof_handler, subdomain_association);
 
-    const unsigned int n_subdomains = 1 + (*std::max_element (subdomain_association.begin (),
-                                                              subdomain_association.end ()   ));
+    const unsigned int n_subdomains
+      = (dynamic_cast<const parallel::Triangulation<DoFHandlerType::dimension,DoFHandlerType::space_dimension> *>
+         (&dof_handler.get_triangulation()) == 0
+         ?
+         1
+         :
+         Utilities::MPI::n_mpi_processes
+         (dynamic_cast<const parallel::Triangulation<DoFHandlerType::dimension,DoFHandlerType::space_dimension> *>
+          (&dof_handler.get_triangulation())->get_communicator()));
+    Assert (n_subdomains >
+            *std::max_element (subdomain_association.begin (),
+                               subdomain_association.end ()),
+            ExcInternalError());
 
-    std::vector<dealii::IndexSet> index_sets (n_subdomains,dealii::IndexSet(dof_handler.n_dofs()));
+    std::vector<dealii::IndexSet> index_sets (n_subdomains,
+                                              dealii::IndexSet(dof_handler.n_dofs()));
 
     // loop over subdomain_association and populate IndexSet when a
     // change in subdomain ID is found
diff --git a/tests/mpi/locally_owned_dofs_per_subdomain.cc b/tests/mpi/locally_owned_dofs_per_subdomain.cc
new file mode 100644 (file)
index 0000000..460119e
--- /dev/null
@@ -0,0 +1,59 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2016 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.
+//
+// ---------------------------------------------------------------------
+
+
+// DoFTools::locally_owned_dofs_per_subdomain wrongly sized the array
+// it returns when we have some processors that don't own any cells
+
+#include "../tests.h"
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/distributed/shared_tria.h>
+#include <deal.II/dofs/dof_handler.h>
+#include <deal.II/fe/fe_q.h>
+
+
+template <int dim>
+void test ()
+{
+  parallel::shared::Triangulation<dim> triangulation(MPI_COMM_WORLD);
+  GridGenerator::hyper_cube(triangulation, -1.0, 1.0);
+
+  FE_Q<dim> fe(1);
+  DoFHandler<dim> dof_handler (triangulation);
+  dof_handler.distribute_dofs(fe);
+
+  // this used to crash here:
+  const IndexSet locally_owned_dofs = dof_handler.locally_owned_dofs();
+
+  deallog << "dim=" << dim
+          << std::endl
+          << "rank=" << Utilities::MPI::this_mpi_process(MPI_COMM_WORLD)
+          << std::endl
+          << "  n_dofs=" << dof_handler.n_dofs()
+          << std::endl
+          << "  n_locally_owned_dofs="
+          << locally_owned_dofs.n_elements()
+          << std::endl;
+}
+
+int main (int argc, char **argv)
+{
+  Utilities::MPI::MPI_InitFinalize mpi_initialization (argc, argv, testing_max_num_threads());
+
+  MPILogInitAll all;
+
+  test<2>();
+  test<3>();
+}
diff --git a/tests/mpi/locally_owned_dofs_per_subdomain.mpirun=2.output b/tests/mpi/locally_owned_dofs_per_subdomain.mpirun=2.output
new file mode 100644 (file)
index 0000000..481ba49
--- /dev/null
@@ -0,0 +1,19 @@
+
+DEAL:0::dim=2
+DEAL:0::rank=0
+DEAL:0::  n_dofs=4
+DEAL:0::  n_locally_owned_dofs=4
+DEAL:0::dim=3
+DEAL:0::rank=0
+DEAL:0::  n_dofs=8
+DEAL:0::  n_locally_owned_dofs=8
+
+DEAL:1::dim=2
+DEAL:1::rank=1
+DEAL:1::  n_dofs=4
+DEAL:1::  n_locally_owned_dofs=0
+DEAL:1::dim=3
+DEAL:1::rank=1
+DEAL:1::  n_dofs=8
+DEAL:1::  n_locally_owned_dofs=0
+

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.