]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
Add tests
authorkronbichler <kronbichler@0785d39b-7218-0410-832d-ea1e28bc413d>
Fri, 31 Jan 2014 20:03:54 +0000 (20:03 +0000)
committerkronbichler <kronbichler@0785d39b-7218-0410-832d-ea1e28bc413d>
Fri, 31 Jan 2014 20:03:54 +0000 (20:03 +0000)
git-svn-id: https://svn.dealii.org/trunk@32361 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/include/deal.II/lac/parallel_vector.templates.h
tests/petsc/copy_to_dealvec.cc
tests/petsc/copy_to_dealvec_block.cc [new file with mode: 0644]
tests/petsc/copy_to_dealvec_block.with_mpi=true.mpirun=4.output [new file with mode: 0644]
tests/trilinos/copy_to_dealvec.cc [new file with mode: 0644]
tests/trilinos/copy_to_dealvec.with_mpi=true.mpirun=4.output [new file with mode: 0644]
tests/trilinos/copy_to_dealvec_block.cc [new file with mode: 0644]
tests/trilinos/copy_to_dealvec_block.with_mpi=true.mpirun=4.output [new file with mode: 0644]

index 60f46ea355e8bea7c10cb612951f599e73ba4f43..252c5b926cadb1abb53e3fe6ca823465e7783fa1 100644 (file)
@@ -248,11 +248,26 @@ namespace parallel
     Vector<Number> &
     Vector<Number>::operator = (const TrilinosWrappers::MPI::Vector &trilinos_vec)
     {
-      Assert(trilinos_vec.locally_owned_elements() == locally_owned_elements(),
-             StandardExceptions::ExcInvalidState());
+      if (trilinos_vec.has_ghost_elements() == false)
+        {
+          Assert(trilinos_vec.locally_owned_elements() == locally_owned_elements(),
+                 StandardExceptions::ExcInvalidState());
+        }
+      else
+        // ghosted trilinos vector must contain the local range of this vector
+        // which is contiguous
+        {
+          Assert((trilinos_vec.locally_owned_elements() & locally_owned_elements())
+                 == locally_owned_elements(),
+                 StandardExceptions::ExcInvalidState());
+        }
 
       // create on trilinos data
-      const VectorView<double> in_view (local_size(), trilinos_vec.begin());
+      const std::size_t start_index =
+        trilinos_vec.vector_partitioner().NumMyElements() > 0 ?
+        trilinos_vec.vector_partitioner().
+        LID(static_cast<TrilinosWrappers::types::int_type>(this->local_range().first)) : 0;
+      const VectorView<double> in_view (local_size(), trilinos_vec.begin()+start_index);
       static_cast<dealii::Vector<Number>&>(vector_view) =
         static_cast<const dealii::Vector<double>&>(in_view);
 
index 2bb2c38f7f7f46aa60f9ee3a34db6da9b3c20c36..08c27f28aa1690b563470d8cf6a7c278773c395b 100644 (file)
@@ -1,7 +1,7 @@
 // ---------------------------------------------------------------------
 // $Id$
 //
-// Copyright (C) 2004 - 2013 by the deal.II authors
+// Copyright (C) 2014 by the deal.II authors
 //
 // This file is part of the deal.II library.
 //
diff --git a/tests/petsc/copy_to_dealvec_block.cc b/tests/petsc/copy_to_dealvec_block.cc
new file mode 100644 (file)
index 0000000..b761d37
--- /dev/null
@@ -0,0 +1,138 @@
+// ---------------------------------------------------------------------
+// $Id$
+//
+// Copyright (C) 2014 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 parallel::distributed::Vector::operator=(PETScWrappers::MPI::BlockVector&)
+
+#include "../tests.h"
+#include <deal.II/lac/petsc_parallel_vector.h>
+#include <deal.II/lac/parallel_vector.h>
+#include <deal.II/lac/parallel_block_vector.h>
+#include <deal.II/lac/petsc_parallel_block_vector.h>
+#include <deal.II/base/index_set.h>
+#include <fstream>
+#include <iostream>
+#include <iomanip>
+#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;
+
+  // each processor owns 2 indices and all
+  // are ghosting Element 1 (the second)
+
+  IndexSet local_active(numproc*2);
+  local_active.add_range(myid*2,myid*2+2);
+  IndexSet local_relevant(numproc*2);
+  local_relevant.add_range(1,2);
+
+  PETScWrappers::MPI::Vector vb_one(local_active, MPI_COMM_WORLD);
+  PETScWrappers::MPI::Vector v_one(local_active, local_relevant, MPI_COMM_WORLD);
+
+  parallel::distributed::Vector<double> copied_one(local_active, local_relevant, MPI_COMM_WORLD);
+
+  // set local values
+  vb_one(myid*2)=myid*2.0;
+  vb_one(myid*2+1)=myid*2.0+1.0;
+
+  vb_one.compress(VectorOperation::insert);
+  vb_one*=2.0;
+  v_one=vb_one;
+
+  PETScWrappers::MPI::BlockVector vb, v;
+  vb.reinit(2);
+  v.reinit(2);
+  parallel::distributed::BlockVector<double> copied(2);
+  for (unsigned int bl=0; bl<2; ++bl)
+    {
+      vb.block(bl) = vb_one;
+      v.block(bl) = v_one;
+      copied.block(bl) = copied_one;
+    }
+  vb.collect_sizes();
+  v.collect_sizes();
+  copied.collect_sizes();
+
+  copied = vb;
+
+  // check local values
+  if (Utilities::MPI::this_mpi_process (MPI_COMM_WORLD) == 0)
+    {
+      deallog << myid*2 << ":" << copied(myid*2) << std::endl;
+      deallog << myid*2+1 << ":" << copied(myid*2+1) << std::endl;
+    }
+
+  for (unsigned int bl=0; bl<2; ++bl)
+    {
+      Assert(copied.block(bl)(myid*2) == myid*4.0, ExcInternalError());
+      Assert(copied.block(bl)(myid*2+1) == myid*4.0+2.0, ExcInternalError());
+    }
+
+  copied = v;
+
+  // check ghost values
+  if (Utilities::MPI::this_mpi_process (MPI_COMM_WORLD) == 0)
+    deallog << "ghost: " << copied(1) << std::endl;
+  Assert(copied(1) == 2.0, ExcInternalError());
+  Assert(copied.block(1)(1) == 2.0, ExcInternalError());
+
+  // check local values
+  if (Utilities::MPI::this_mpi_process (MPI_COMM_WORLD) == 0)
+    {
+      deallog << myid*2 << ":" << copied(myid*2) << std::endl;
+      deallog << myid*2+1 << ":" << copied(myid*2+1) << std::endl;
+    }
+
+  for (unsigned int bl=0; bl<2; ++bl)
+    {
+      Assert(copied.block(bl)(myid*2) == myid*4.0, ExcInternalError());
+      Assert(copied.block(bl)(myid*2+1) == myid*4.0+2.0, ExcInternalError());
+    }
+
+  // done
+  if (myid==0)
+    deallog << "OK" << std::endl;
+}
+
+
+
+int main (int argc, char **argv)
+{
+  Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv);
+  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/petsc/copy_to_dealvec_block.with_mpi=true.mpirun=4.output b/tests/petsc/copy_to_dealvec_block.with_mpi=true.mpirun=4.output
new file mode 100644 (file)
index 0000000..6ed3f01
--- /dev/null
@@ -0,0 +1,8 @@
+
+DEAL:0::numproc=4
+DEAL:0::0:0
+DEAL:0::1:2.000
+DEAL:0::ghost: 2.000
+DEAL:0::0:0
+DEAL:0::1:2.000
+DEAL:0::OK
diff --git a/tests/trilinos/copy_to_dealvec.cc b/tests/trilinos/copy_to_dealvec.cc
new file mode 100644 (file)
index 0000000..a723463
--- /dev/null
@@ -0,0 +1,120 @@
+// ---------------------------------------------------------------------
+// $Id$
+//
+// Copyright (C) 2014 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 parallel::distributed::Vector::operator=(TrilinosWrappers::MPI::Vector&)
+
+#include "../tests.h"
+#include <deal.II/lac/trilinos_vector.h>
+#include <deal.II/lac/parallel_vector.h>
+#include <deal.II/base/index_set.h>
+#include <fstream>
+#include <iostream>
+#include <iomanip>
+#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;
+
+  // each processor owns 2 indices and all
+  // are ghosting Element 1 (the second)
+
+  IndexSet local_active(numproc*2);
+  local_active.add_range(myid*2,myid*2+2);
+  IndexSet local_relevant(numproc*2);
+  local_relevant.add_range(1,2);
+
+  TrilinosWrappers::MPI::Vector vb(local_active, MPI_COMM_WORLD);
+  TrilinosWrappers::MPI::Vector v(local_active, local_relevant, MPI_COMM_WORLD);
+
+  parallel::distributed::Vector<double> copied(local_active, local_relevant, MPI_COMM_WORLD);
+
+  // set local values
+  vb(myid*2)=myid*2.0;
+  vb(myid*2+1)=myid*2.0+1.0;
+
+  vb.compress(VectorOperation::insert);
+  vb*=2.0;
+  v=vb;
+  //v.update_ghost_values();
+
+  Assert(!vb.has_ghost_elements(), ExcInternalError());
+  Assert(v.has_ghost_elements(), ExcInternalError());
+
+  copied = vb;
+
+  // check local values
+  if (Utilities::MPI::this_mpi_process (MPI_COMM_WORLD) == 0)
+    {
+      deallog << myid*2 << ":" << copied(myid*2) << std::endl;
+      deallog << myid*2+1 << ":" << copied(myid*2+1) << std::endl;
+    }
+
+  Assert(copied(myid*2) == myid*4.0, ExcInternalError());
+  Assert(copied(myid*2+1) == myid*4.0+2.0, ExcInternalError());
+
+  copied = v;
+
+  // check ghost values
+  if (Utilities::MPI::this_mpi_process (MPI_COMM_WORLD) == 0)
+    deallog << "ghost: " << copied(1) << std::endl;
+  Assert(copied(1) == 2.0, ExcInternalError());
+
+  // check local values
+  if (Utilities::MPI::this_mpi_process (MPI_COMM_WORLD) == 0)
+    {
+      deallog << myid*2 << ":" << copied(myid*2) << std::endl;
+      deallog << myid*2+1 << ":" << copied(myid*2+1) << std::endl;
+    }
+
+  Assert(copied(myid*2) == myid*4.0, ExcInternalError());
+  Assert(copied(myid*2+1) == myid*4.0+2.0, ExcInternalError());
+
+
+  // done
+  if (myid==0)
+    deallog << "OK" << std::endl;
+}
+
+
+
+int main (int argc, char **argv)
+{
+  Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv);
+  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/trilinos/copy_to_dealvec.with_mpi=true.mpirun=4.output b/tests/trilinos/copy_to_dealvec.with_mpi=true.mpirun=4.output
new file mode 100644 (file)
index 0000000..6ed3f01
--- /dev/null
@@ -0,0 +1,8 @@
+
+DEAL:0::numproc=4
+DEAL:0::0:0
+DEAL:0::1:2.000
+DEAL:0::ghost: 2.000
+DEAL:0::0:0
+DEAL:0::1:2.000
+DEAL:0::OK
diff --git a/tests/trilinos/copy_to_dealvec_block.cc b/tests/trilinos/copy_to_dealvec_block.cc
new file mode 100644 (file)
index 0000000..15b41fe
--- /dev/null
@@ -0,0 +1,136 @@
+// ---------------------------------------------------------------------
+// $Id$
+//
+// Copyright (C) 2014 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 parallel::distributed::Vector::operator=(TrilinosWrappers::MPI::BlockVector&)
+
+#include "../tests.h"
+#include <deal.II/lac/trilinos_vector.h>
+#include <deal.II/lac/parallel_vector.h>
+#include <deal.II/lac/parallel_block_vector.h>
+#include <deal.II/lac/trilinos_parallel_block_vector.h>
+#include <deal.II/base/index_set.h>
+#include <fstream>
+#include <iostream>
+#include <iomanip>
+#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;
+
+  // each processor owns 2 indices and all
+  // are ghosting Element 1 (the second)
+
+  IndexSet local_active(numproc*2);
+  local_active.add_range(myid*2,myid*2+2);
+  IndexSet local_relevant(numproc*2);
+  local_relevant.add_range(1,2);
+
+  TrilinosWrappers::MPI::Vector vb_one(local_active, MPI_COMM_WORLD);
+  TrilinosWrappers::MPI::Vector v_one(local_active, local_relevant, MPI_COMM_WORLD);
+
+  parallel::distributed::Vector<double> copied_one(local_active, local_relevant, MPI_COMM_WORLD);
+
+  // set local values
+  vb_one(myid*2)=myid*2.0;
+  vb_one(myid*2+1)=myid*2.0+1.0;
+
+  vb_one.compress(VectorOperation::insert);
+  vb_one*=2.0;
+  v_one=vb_one;
+
+  TrilinosWrappers::MPI::BlockVector vb(2), v(2);
+  parallel::distributed::BlockVector<double> copied(2);
+  for (unsigned int bl=0; bl<2; ++bl)
+    {
+      vb.block(bl) = vb_one;
+      v.block(bl) = v_one;
+      copied.block(bl) = copied_one;
+    }
+  vb.collect_sizes();
+  v.collect_sizes();
+  copied.collect_sizes();
+
+  copied = vb;
+
+  // check local values
+  if (Utilities::MPI::this_mpi_process (MPI_COMM_WORLD) == 0)
+    {
+      deallog << myid*2 << ":" << copied(myid*2) << std::endl;
+      deallog << myid*2+1 << ":" << copied(myid*2+1) << std::endl;
+    }
+
+  for (unsigned int bl=0; bl<2; ++bl)
+    {
+      Assert(copied.block(bl)(myid*2) == myid*4.0, ExcInternalError());
+      Assert(copied.block(bl)(myid*2+1) == myid*4.0+2.0, ExcInternalError());
+    }
+
+  copied = v;
+
+  // check ghost values
+  if (Utilities::MPI::this_mpi_process (MPI_COMM_WORLD) == 0)
+    deallog << "ghost: " << copied(1) << std::endl;
+  Assert(copied(1) == 2.0, ExcInternalError());
+  Assert(copied.block(1)(1) == 2.0, ExcInternalError());
+
+  // check local values
+  if (Utilities::MPI::this_mpi_process (MPI_COMM_WORLD) == 0)
+    {
+      deallog << myid*2 << ":" << copied(myid*2) << std::endl;
+      deallog << myid*2+1 << ":" << copied(myid*2+1) << std::endl;
+    }
+
+  for (unsigned int bl=0; bl<2; ++bl)
+    {
+      Assert(copied.block(bl)(myid*2) == myid*4.0, ExcInternalError());
+      Assert(copied.block(bl)(myid*2+1) == myid*4.0+2.0, ExcInternalError());
+    }
+
+  // done
+  if (myid==0)
+    deallog << "OK" << std::endl;
+}
+
+
+
+int main (int argc, char **argv)
+{
+  Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv);
+  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/trilinos/copy_to_dealvec_block.with_mpi=true.mpirun=4.output b/tests/trilinos/copy_to_dealvec_block.with_mpi=true.mpirun=4.output
new file mode 100644 (file)
index 0000000..6ed3f01
--- /dev/null
@@ -0,0 +1,8 @@
+
+DEAL:0::numproc=4
+DEAL:0::0:0
+DEAL:0::1:2.000
+DEAL:0::ghost: 2.000
+DEAL:0::0:0
+DEAL:0::1:2.000
+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.