]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Fix a bug in make_trilinos_map() where different processors would call different... 2242/head
authorBruno Turcksin <bruno.turcksin@gmail.com>
Wed, 24 Feb 2016 20:10:59 +0000 (15:10 -0500)
committerBruno Turcksin <bruno.turcksin@gmail.com>
Wed, 24 Feb 2016 22:57:36 +0000 (17:57 -0500)
doc/news/changes.h
source/base/index_set.cc
tests/trilinos/index_set_01.cc [new file with mode: 0644]
tests/trilinos/index_set_01.mpirun=3.output [new file with mode: 0644]

index 1f84ee4d2df5be52f87b9c6bbd0a4e2fb066ec25..7d2f1790a2eac179e06b2b5413f1852e6eeb3a84 100644 (file)
@@ -63,6 +63,12 @@ inconvenience this causes.
 
 
 <ol> 
+ <li> Fixed: The function IndexSet::make_trilinos_map() now works if some 
+ processors have a contiguous range of indices and others do not.
+ <br>
+ (Bruno Turcksin, 2016/02/17)
+ </li>
+
  <li> Updated: step-44 has been been expressed in a more dimension independent 
  manner, and can be now run in both 2-d and 3-d.
  <br>
index 3cc560614b93fc62217b5997064dcc21363f366e..5850d56eb0f7b8b8780dd5889c6052289a27b54d 100644 (file)
@@ -510,7 +510,17 @@ IndexSet::make_trilinos_map (const MPI_Comm &communicator,
     }
 #endif
 
+  // Check that all the processors have a contiguous range of values. Otherwise,
+  // we risk to call different Epetra_Map on different processors and the code
+  // hangs.
+#ifdef DEAL_II_WITH_MPI
+  bool all_contiguous = is_contiguous();
+  MPI_Allreduce(MPI_IN_PLACE, &all_contiguous, 1, MPI_LOGICAL, MPI_LAND,
+                communicator);
+  if ((all_contiguous == true) && (!overlapping))
+#else
   if ((is_contiguous() == true) && (!overlapping))
+#endif
     return Epetra_Map (TrilinosWrappers::types::int_type(size()),
                        TrilinosWrappers::types::int_type(n_elements()),
                        0,
diff --git a/tests/trilinos/index_set_01.cc b/tests/trilinos/index_set_01.cc
new file mode 100644 (file)
index 0000000..574f399
--- /dev/null
@@ -0,0 +1,66 @@
+// ---------------------------------------------------------------------
+//
+// 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.
+//
+// ---------------------------------------------------------------------
+
+
+
+// check that make_trilinos_map always uses the same constructor all the
+// processor. If we don't the code will hang like reported on the mailing list
+// by Nicola Giuliani.
+
+#include "../tests.h"
+#include <deal.II/base/utilities.h>
+#include <deal.II/base/index_set.h>
+#include <fstream>
+#include <iostream>
+
+
+
+void test ()
+{
+  IndexSet my_set;
+  my_set.set_size(654);
+  int this_mpi_process, n_mpi_processes;
+  n_mpi_processes = Utilities::MPI::n_mpi_processes(MPI_COMM_WORLD);
+  this_mpi_process = Utilities::MPI::this_mpi_process(MPI_COMM_WORLD);
+  if (this_mpi_process == 0)
+    {
+      my_set.add_range(0,86);
+    }
+  else if (this_mpi_process == 1)
+    {
+      my_set.add_range(86,400);
+      my_set.add_range(529,654);
+    }
+  else
+    {
+      my_set.add_range(400,529);
+    }
+  my_set.compress();
+  my_set.make_trilinos_map(MPI_COMM_WORLD,false);
+}
+
+
+
+int main (int argc,char **argv)
+{
+  std::ofstream logfile("output");
+  deallog.attach(logfile);
+  deallog.threshold_double(1.e-10);
+
+  Utilities::MPI::MPI_InitFinalize mpi_initialization (argc, argv, testing_max_num_threads());
+  test();
+
+  deallog << "OK" << std::endl;
+}
diff --git a/tests/trilinos/index_set_01.mpirun=3.output b/tests/trilinos/index_set_01.mpirun=3.output
new file mode 100644 (file)
index 0000000..0fd8fc1
--- /dev/null
@@ -0,0 +1,2 @@
+
+DEAL::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.