]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Fix a bug in import function of distributed::Vector. 3859/head
authorBruno Turcksin <bruno.turcksin@gmail.com>
Fri, 27 Jan 2017 14:50:32 +0000 (09:50 -0500)
committerBruno Turcksin <bruno.turcksin@gmail.com>
Fri, 27 Jan 2017 14:50:32 +0000 (09:50 -0500)
include/deal.II/lac/la_parallel_vector.templates.h
include/deal.II/lac/read_write_vector.templates.h
tests/mpi/parallel_vector_20.cc [new file with mode: 0644]
tests/mpi/parallel_vector_20.mpirun=2.output [new file with mode: 0644]

index 4b8fdd1b5114143b33436c3a523a0d24aabab289..2033688deaa737820edb0b96ca99efd265fa33b9 100644 (file)
@@ -827,13 +827,20 @@ namespace LinearAlgebra
                            VectorOperation::values                         operation,
                            std_cxx11::shared_ptr<const CommunicationPatternBase> communication_pattern)
     {
+      IndexSet locally_owned_elem = locally_owned_elements();
       // If no communication pattern is given, create one. Otherwise, use the
       // given one.
       std_cxx11::shared_ptr<const Utilities::MPI::Partitioner> comm_pattern;
       if (communication_pattern.get() == NULL)
         {
-          comm_pattern.reset(new Utilities::MPI::Partitioner(locally_owned_elements(),
-                                                             V.get_stored_elements(),
+          // Split the IndexSet of V in locally owned elements and ghost indices
+          // then create the communication pattern
+          IndexSet ghost_indices(V.get_stored_elements());
+          ghost_indices.subtract_set(locally_owned_elem);
+          IndexSet local_indices(V.get_stored_elements());
+          local_indices.subtract_set(ghost_indices);
+          comm_pattern.reset(new Utilities::MPI::Partitioner(local_indices,
+                                                             ghost_indices,
                                                              get_mpi_communicator()));
         }
       else
@@ -856,9 +863,13 @@ namespace LinearAlgebra
 
       tmp_vector.compress(operation);
 
-      dealii::internal::VectorOperations::Vector_copy<Number,Number> copier(tmp_vector.val, val);
-      internal::VectorOperations::parallel_for(copier, partitioner->local_size(),
-                                               thread_loop_partitioner);
+      // Copy the local elements of tmp_vector to the right place in val
+      IndexSet tmp_index_set = tmp_vector.locally_owned_elements();
+      for (size_type i=0; i<tmp_index_set.n_elements(); ++i)
+        {
+          val[locally_owned_elem.index_within_set(tmp_index_set.nth_index_in_set(i))] =
+            tmp_vector.local_element(i);
+        }
     }
 
 
index f24c76ced18a5360308aaec4a706c84acb85074f..2f46df57ebfa3e50ed13f8c8a01a7d20056a0d9c 100644 (file)
@@ -126,7 +126,7 @@ namespace LinearAlgebra
     if (omit_zeroing_entries == false)
       this->operator= (Number());
 
-    // reset the communication patter
+    // reset the communication pattern
     source_stored_elements.clear();
     comm_pattern.reset();
   }
diff --git a/tests/mpi/parallel_vector_20.cc b/tests/mpi/parallel_vector_20.cc
new file mode 100644 (file)
index 0000000..fb9ed93
--- /dev/null
@@ -0,0 +1,84 @@
+// ---------------------------------------------------------------------
+//
+// 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.
+//
+// ---------------------------------------------------------------------
+
+// Check import function
+
+#include "../tests.h"
+#include <deal.II/base/utilities.h>
+#include <deal.II/base/index_set.h>
+#include <deal.II/lac/la_parallel_vector.h>
+#include <deal.II/lac/read_write_vector.h>
+#include <fstream>
+#include <iostream>
+#include <vector>
+
+void test()
+{
+  unsigned int my_id = Utilities::MPI::this_mpi_process (MPI_COMM_WORLD);
+  unsigned int n_procs = Utilities::MPI::n_mpi_processes (MPI_COMM_WORLD);
+
+  IndexSet locally_owned(n_procs*2);
+  locally_owned.add_range(my_id*2, my_id*2+2);
+  IndexSet read_write_owned(4);
+  read_write_owned.add_range(my_id*2, my_id*2+2);
+
+  LinearAlgebra::distributed::Vector<double> v(locally_owned, MPI_COMM_WORLD);
+  LinearAlgebra::ReadWriteVector<double> read_write_vector(read_write_owned);
+  read_write_vector.local_element(0) = 1.;
+  read_write_vector.local_element(1) = 2.;
+
+  v.import(read_write_vector, VectorOperation::values::insert);
+
+  AssertThrow(v.local_element(0) == 1., ExcInternalError());
+  AssertThrow(v.local_element(1) == 2., ExcInternalError());
+
+  read_write_owned.clear();
+  read_write_owned.add_index(1);
+  read_write_owned.add_index(2);
+  read_write_vector.reinit(read_write_owned);
+  read_write_vector.local_element(0) = 1.;
+  read_write_vector.local_element(1) = 2.;
+
+  v.import(read_write_vector, VectorOperation::values::insert);
+
+  AssertThrow(v.local_element(0) == my_id+1, ExcInternalError());
+  AssertThrow(v.local_element(1) == my_id+1, ExcInternalError());
+
+  if (my_id == 0)
+    deallog << "OK" << std::endl;
+}
+
+
+
+int main (int argc, char **argv)
+{
+  Utilities::MPI::MPI_InitFinalize mpi_initialization (argc, argv, testing_max_num_threads());
+
+  unsigned int myid = Utilities::MPI::this_mpi_process (MPI_COMM_WORLD);
+  deallog.push(Utilities::int_to_string(myid));
+
+  if (myid == 0)
+    {
+      std::ofstream logfile("output");
+      deallog.attach(logfile);
+      deallog << std::setprecision(4);
+      deallog.threshold_double(1.e-10);
+
+      test();
+    }
+  else
+    test();
+
+}
diff --git a/tests/mpi/parallel_vector_20.mpirun=2.output b/tests/mpi/parallel_vector_20.mpirun=2.output
new file mode 100644 (file)
index 0000000..be8d055
--- /dev/null
@@ -0,0 +1,2 @@
+
+DEAL:0::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.