]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Bugfix in state 'has_ghost_elements' of parallel vector. 419/head
authorMartin Kronbichler <kronbichler@lnm.mw.tum.de>
Mon, 12 Jan 2015 09:13:41 +0000 (10:13 +0100)
committerMartin Kronbichler <kronbichler@lnm.mw.tum.de>
Mon, 12 Jan 2015 09:17:20 +0000 (10:17 +0100)
When assigning a non-ghosted vector from another non-ghosted vector,
the state of the vector should be non-ghosted also when the partitioners
of the vectors are different.

include/deal.II/lac/parallel_vector.h
tests/mpi/parallel_vector_15.cc [new file with mode: 0644]
tests/mpi/parallel_vector_15.mpirun=1.output [new file with mode: 0644]
tests/mpi/parallel_vector_15.mpirun=3.output [new file with mode: 0644]

index 457782cfb8d8101fc5d0a12a6e5053c1d7b17c65..ad8e2d5f7cf54a2186fb85113fae195ceed81dc9 100644 (file)
@@ -1,6 +1,6 @@
 // ---------------------------------------------------------------------
 //
-// Copyright (C) 2011 - 2014 by the deal.II authors
+// Copyright (C) 2011 - 2015 by the deal.II authors
 //
 // This file is part of the deal.II library.
 //
@@ -1252,7 +1252,7 @@ namespace parallel
       // we update ghost values whenever one of the input or output vector
       // already held ghost values or when we import data from a vector with
       // the same local range but different ghost layout
-      bool must_update_ghost_values = true;
+      bool must_update_ghost_values = c.vector_is_ghosted;
 
       // check whether the two vectors use the same parallel partitioner. if
       // not, check if all local ranges are the same (that way, we can
@@ -1269,9 +1269,11 @@ namespace parallel
               ||
               local_ranges_different_loc)
             reinit (c, true);
+          else
+            must_update_ghost_values |= vector_is_ghosted;
         }
       else
-        must_update_ghost_values = vector_is_ghosted || c.vector_is_ghosted;
+        must_update_ghost_values |= vector_is_ghosted;
 
       vector_view = c.vector_view;
       if (must_update_ghost_values)
diff --git a/tests/mpi/parallel_vector_15.cc b/tests/mpi/parallel_vector_15.cc
new file mode 100644 (file)
index 0000000..7117d38
--- /dev/null
@@ -0,0 +1,117 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2011 - 2015 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 handling of ghost elements in parallel distributed vectors works
+// appropriately when creating a vector from a non-ghosted source vector using
+// the assignment operator
+
+#include "../tests.h"
+#include <deal.II/base/utilities.h>
+#include <deal.II/base/index_set.h>
+#include <deal.II/lac/parallel_vector.h>
+#include <fstream>
+#include <iostream>
+#include <vector>
+
+
+void test ()
+{
+  unsigned int myid = Utilities::MPI::this_mpi_process (MPI_COMM_WORLD);
+  unsigned int numproc = Utilities::MPI::n_mpi_processes (MPI_COMM_WORLD);
+
+  if (myid==0) deallog << "numproc=" << numproc << std::endl;
+
+  // processor 0 and 1 own 2 indices each, higher processors nothing, all are
+  // ghosting global elements 1 and 3
+  IndexSet local_owned(std::min(numproc*2, 4U));
+  if (myid < 2)
+    local_owned.add_range(myid*2,myid*2+2);
+  IndexSet local_relevant(local_owned.size());
+  local_relevant = local_owned;
+  local_relevant.add_range(1,2);
+  if (numproc > 1)
+    local_relevant.add_range(3,4);
+
+  parallel::distributed::Vector<double> v(local_owned, local_relevant, MPI_COMM_WORLD);
+
+  // set local values
+  if (myid < 2)
+    {
+      v(myid*2)=myid*2.0;
+      v(myid*2+1)=myid*2.0+1.0;
+    }
+
+  v.compress(VectorOperation::insert);
+
+  if (myid==0)
+    deallog << "v has ghost elements: " << v.has_ghost_elements() << std::endl;
+
+  parallel::distributed::Vector<double> w, x;
+  w = v;
+  if (myid==0)
+    deallog << "w has ghost elements: " << w.has_ghost_elements() << std::endl;
+
+  v.update_ghost_values();
+  w = v;
+  if (myid==0)
+    deallog << "w has ghost elements: " << w.has_ghost_elements() << std::endl;
+
+  v.zero_out_ghosts();
+  w = v;
+  if (myid==0)
+    deallog << "w has ghost elements: " << w.has_ghost_elements() << std::endl;
+
+  w.zero_out_ghosts();
+  w = v;
+  if (myid==0)
+    deallog << "w has ghost elements: " << w.has_ghost_elements() << std::endl;
+
+  v.update_ghost_values();
+  x = v;
+  if (myid==0)
+    deallog << "x has ghost elements: " << x.has_ghost_elements() << std::endl;
+
+  x.zero_out_ghosts();
+  if (myid==0)
+    deallog << "x has ghost elements: " << x.has_ghost_elements() << std::endl;
+
+  if (myid==0)
+    deallog << "OK" << std::endl;
+}
+
+
+
+int main (int argc, char **argv)
+{
+  Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, numbers::invalid_unsigned_int);
+
+  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.depth_console(0);
+      deallog.threshold_double(1.e-10);
+
+      test();
+    }
+  else
+    test();
+
+}
diff --git a/tests/mpi/parallel_vector_15.mpirun=1.output b/tests/mpi/parallel_vector_15.mpirun=1.output
new file mode 100644 (file)
index 0000000..a70fb56
--- /dev/null
@@ -0,0 +1,10 @@
+
+DEAL:0::numproc=1
+DEAL:0::v has ghost elements: 0
+DEAL:0::w has ghost elements: 0
+DEAL:0::w has ghost elements: 1
+DEAL:0::w has ghost elements: 1
+DEAL:0::w has ghost elements: 0
+DEAL:0::x has ghost elements: 1
+DEAL:0::x has ghost elements: 0
+DEAL:0::OK
diff --git a/tests/mpi/parallel_vector_15.mpirun=3.output b/tests/mpi/parallel_vector_15.mpirun=3.output
new file mode 100644 (file)
index 0000000..2588b1d
--- /dev/null
@@ -0,0 +1,10 @@
+
+DEAL:0::numproc=3
+DEAL:0::v has ghost elements: 0
+DEAL:0::w has ghost elements: 0
+DEAL:0::w has ghost elements: 1
+DEAL:0::w has ghost elements: 1
+DEAL:0::w has ghost elements: 0
+DEAL:0::x has ghost elements: 1
+DEAL:0::x has ghost elements: 0
+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.