]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add tests for native deal variants and TrilinosWrappers::MPI
authorMatthias Maier <tamiko@43-1.org>
Fri, 1 May 2015 22:50:33 +0000 (00:50 +0200)
committerMatthias Maier <tamiko@43-1.org>
Sat, 2 May 2015 10:18:11 +0000 (12:18 +0200)
tests/lac/vector_move_operations.cc [new file with mode: 0644]
tests/lac/vector_move_operations.with_cxx11=on.output [new file with mode: 0644]
tests/trilinos/vector_move_operations.cc [new file with mode: 0644]
tests/trilinos/vector_move_operations.with_cxx11=on.output [new file with mode: 0644]

diff --git a/tests/lac/vector_move_operations.cc b/tests/lac/vector_move_operations.cc
new file mode 100644 (file)
index 0000000..9fbc513
--- /dev/null
@@ -0,0 +1,92 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2012 - 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.
+//
+// ---------------------------------------------------------------------
+
+
+#include "../tests.h"
+#include <deal.II/lac/vector.h>
+#include <deal.II/lac/block_vector.h>
+
+#define PRINTME(name, var) \
+  deallog << "Block vector: " name << ":" << std::endl; \
+  for (unsigned int i = 0; i < var.n_blocks(); ++i) \
+    deallog << "[block " << i << " ]  " << var.block(i);
+
+int main()
+{
+  initlog();
+
+  {
+    Vector<double> w(Vector<double>(10));
+    deallog << "move constructor:  " << w;
+  }
+
+  deallog << std::endl;
+
+  {
+    Vector<double> u(10);
+    for (unsigned int i = 0; i < u.size(); ++i)
+      u[i] = (double)(i+1);
+
+    deallog << "vector:          " << u;
+
+    Vector<double> v;
+    v = u;
+    deallog << "copy assignment: " << v;
+    deallog << "old object:      " << u;
+
+    v = std::move(u);
+    deallog << "move assignment: " << v;
+    deallog << "old object size: " << u.size() << std::endl;
+
+    // and swap again with different sizes
+    u = std::move(v);
+    deallog << "move assignment: " << u;
+    deallog << "old object size: " << v.size() << std::endl;
+  }
+
+  deallog << std::endl;
+
+  {
+    BlockVector<double> w(BlockVector<double>(5, 2));
+    PRINTME("move constructor", w);
+  }
+
+  deallog << std::endl;
+
+  {
+    BlockVector<double> u(5, 2);
+    for (unsigned int i = 0; i < 5; ++i)
+      for (unsigned int j = 0; j < 2; ++j)
+        u.block(i)[j] = (double)(10 * i + j);
+
+    PRINTME("BlockVector", u);
+
+    BlockVector<double> v;
+    v = u;
+    PRINTME("copy assignemnt", v);
+    PRINTME("old object", u);
+
+    v = std::move(u);
+    PRINTME("move assignemnt", v);
+    deallog << "old object size: " << u.n_blocks() << std::endl;
+
+    // and swap again with different sizes
+    u = std::move(v);
+    PRINTME("move assignemnt", u);
+    deallog << "old object size: " << v.n_blocks() << std::endl;
+  }
+}
+
+
diff --git a/tests/lac/vector_move_operations.with_cxx11=on.output b/tests/lac/vector_move_operations.with_cxx11=on.output
new file mode 100644 (file)
index 0000000..29f0014
--- /dev/null
@@ -0,0 +1,50 @@
+
+DEAL::move constructor:  0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 
+DEAL::
+DEAL::vector:          1.00000 2.00000 3.00000 4.00000 5.00000 6.00000 7.00000 8.00000 9.00000 10.0000 
+DEAL::copy assignment: 1.00000 2.00000 3.00000 4.00000 5.00000 6.00000 7.00000 8.00000 9.00000 10.0000 
+DEAL::old object:      1.00000 2.00000 3.00000 4.00000 5.00000 6.00000 7.00000 8.00000 9.00000 10.0000 
+DEAL::move assignment: 1.00000 2.00000 3.00000 4.00000 5.00000 6.00000 7.00000 8.00000 9.00000 10.0000 
+DEAL::old object size: 0
+DEAL::move assignment: 1.00000 2.00000 3.00000 4.00000 5.00000 6.00000 7.00000 8.00000 9.00000 10.0000 
+DEAL::old object size: 0
+DEAL::
+DEAL::Block vector: move constructor:
+DEAL::[block 0 ]  0.00000 0.00000 
+DEAL::[block 1 ]  0.00000 0.00000 
+DEAL::[block 2 ]  0.00000 0.00000 
+DEAL::[block 3 ]  0.00000 0.00000 
+DEAL::[block 4 ]  0.00000 0.00000 
+DEAL::
+DEAL::Block vector: BlockVector:
+DEAL::[block 0 ]  0.00000 1.00000 
+DEAL::[block 1 ]  10.0000 11.0000 
+DEAL::[block 2 ]  20.0000 21.0000 
+DEAL::[block 3 ]  30.0000 31.0000 
+DEAL::[block 4 ]  40.0000 41.0000 
+DEAL::Block vector: copy assignemnt:
+DEAL::[block 0 ]  0.00000 1.00000 
+DEAL::[block 1 ]  10.0000 11.0000 
+DEAL::[block 2 ]  20.0000 21.0000 
+DEAL::[block 3 ]  30.0000 31.0000 
+DEAL::[block 4 ]  40.0000 41.0000 
+DEAL::Block vector: old object:
+DEAL::[block 0 ]  0.00000 1.00000 
+DEAL::[block 1 ]  10.0000 11.0000 
+DEAL::[block 2 ]  20.0000 21.0000 
+DEAL::[block 3 ]  30.0000 31.0000 
+DEAL::[block 4 ]  40.0000 41.0000 
+DEAL::Block vector: move assignemnt:
+DEAL::[block 0 ]  0.00000 1.00000 
+DEAL::[block 1 ]  10.0000 11.0000 
+DEAL::[block 2 ]  20.0000 21.0000 
+DEAL::[block 3 ]  30.0000 31.0000 
+DEAL::[block 4 ]  40.0000 41.0000 
+DEAL::old object size: 0
+DEAL::Block vector: move assignemnt:
+DEAL::[block 0 ]  0.00000 1.00000 
+DEAL::[block 1 ]  10.0000 11.0000 
+DEAL::[block 2 ]  20.0000 21.0000 
+DEAL::[block 3 ]  30.0000 31.0000 
+DEAL::[block 4 ]  40.0000 41.0000 
+DEAL::old object size: 0
diff --git a/tests/trilinos/vector_move_operations.cc b/tests/trilinos/vector_move_operations.cc
new file mode 100644 (file)
index 0000000..1a57fba
--- /dev/null
@@ -0,0 +1,92 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2012 - 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.
+//
+// ---------------------------------------------------------------------
+
+
+#include "../tests.h"
+#include <deal.II/lac/trilinos_vector.h>
+#include <deal.II/lac/trilinos_parallel_block_vector.h>
+
+#define PRINTME(name, var)                                 \
+  deallog << "Vector: " name << ":" << std::endl;          \
+  for (unsigned int i = 0; i < var.size(); ++i)            \
+    deallog << var[i] << " ";                              \
+  deallog << std::endl;
+
+#define PRINTBLOCK(name, var)                                \
+  deallog << "Block vector: " name << ":" << std::endl;      \
+  for (unsigned int i = 0; i < var.n_blocks(); ++i)          \
+    {                                                        \
+      deallog << "[block " << i << " ]  ";                   \
+      for (unsigned int j = 0; j < var.block(i).size(); ++j) \
+        deallog << var.block(i)[j] << " ";                   \
+      deallog << std::endl;                                  \
+    }
+
+
+int main (int argc, char **argv)
+{
+  initlog();
+
+  Utilities::MPI::MPI_InitFinalize mpi_initialization (argc, argv, numbers::invalid_unsigned_int);
+
+  {
+    IndexSet local_owned(10);
+    local_owned.add_range (0,10);
+
+    TrilinosWrappers::MPI::Vector u(local_owned, MPI_COMM_WORLD);
+    for (unsigned int i = 0; i < u.size(); ++i)
+      u[i] = (double)(i+1);
+
+    PRINTME("Vector", u);
+
+    TrilinosWrappers::MPI::Vector v;
+    v = u;
+    PRINTME("copy assignemnt", v);
+    PRINTME("old object", u);
+
+    v = std::move(u);
+    PRINTME("move assignemnt", v);
+    deallog << "old object size: " << u.size() << std::endl;
+  }
+
+  deallog << std::endl;
+
+  {
+    std::vector<IndexSet> local_owned(5);
+    for (auto &index : local_owned)
+      {
+        index.set_size(2);
+        index.add_range(0,2);
+      }
+
+    TrilinosWrappers::MPI::BlockVector u(local_owned, MPI_COMM_WORLD);
+    for (unsigned int i = 0; i < 5; ++i)
+      for (unsigned int j = 0; j < 2; ++j)
+        u.block(i)[j] = (double)(10 * i + j);
+
+    PRINTBLOCK("BlockVector", u);
+
+    TrilinosWrappers::MPI::BlockVector v;
+    v = u;
+    PRINTBLOCK("copy assignemnt", v);
+    PRINTBLOCK("old object", u);
+
+    v = std::move(u);
+    PRINTBLOCK("move assignemnt", v);
+    deallog << "old object size: " << u.n_blocks() << std::endl;
+  }
+}
+
+
diff --git a/tests/trilinos/vector_move_operations.with_cxx11=on.output b/tests/trilinos/vector_move_operations.with_cxx11=on.output
new file mode 100644 (file)
index 0000000..ab9bb97
--- /dev/null
@@ -0,0 +1,36 @@
+
+DEAL::Vector: Vector:
+DEAL::1.00000 2.00000 3.00000 4.00000 5.00000 6.00000 7.00000 8.00000 9.00000 10.0000 
+DEAL::Vector: copy assignemnt:
+DEAL::1.00000 2.00000 3.00000 4.00000 5.00000 6.00000 7.00000 8.00000 9.00000 10.0000 
+DEAL::Vector: old object:
+DEAL::1.00000 2.00000 3.00000 4.00000 5.00000 6.00000 7.00000 8.00000 9.00000 10.0000 
+DEAL::Vector: move assignemnt:
+DEAL::1.00000 2.00000 3.00000 4.00000 5.00000 6.00000 7.00000 8.00000 9.00000 10.0000 
+DEAL::old object size: 0
+DEAL::
+DEAL::Block vector: BlockVector:
+DEAL::[block 0 ]  0.00000 1.00000 
+DEAL::[block 1 ]  10.0000 11.0000 
+DEAL::[block 2 ]  20.0000 21.0000 
+DEAL::[block 3 ]  30.0000 31.0000 
+DEAL::[block 4 ]  40.0000 41.0000 
+DEAL::Block vector: copy assignemnt:
+DEAL::[block 0 ]  0.00000 1.00000 
+DEAL::[block 1 ]  10.0000 11.0000 
+DEAL::[block 2 ]  20.0000 21.0000 
+DEAL::[block 3 ]  30.0000 31.0000 
+DEAL::[block 4 ]  40.0000 41.0000 
+DEAL::Block vector: old object:
+DEAL::[block 0 ]  0.00000 1.00000 
+DEAL::[block 1 ]  10.0000 11.0000 
+DEAL::[block 2 ]  20.0000 21.0000 
+DEAL::[block 3 ]  30.0000 31.0000 
+DEAL::[block 4 ]  40.0000 41.0000 
+DEAL::Block vector: move assignemnt:
+DEAL::[block 0 ]  0.00000 1.00000 
+DEAL::[block 1 ]  10.0000 11.0000 
+DEAL::[block 2 ]  20.0000 21.0000 
+DEAL::[block 3 ]  30.0000 31.0000 
+DEAL::[block 4 ]  40.0000 41.0000 
+DEAL::old object size: 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.