]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Added move constructor to SparseMatrix
authordanshapero <shapero.daniel@gmail.com>
Thu, 24 Mar 2016 00:47:10 +0000 (17:47 -0700)
committerdanshapero <shapero.daniel@gmail.com>
Thu, 24 Mar 2016 00:47:10 +0000 (17:47 -0700)
include/deal.II/lac/sparse_matrix.h
include/deal.II/lac/sparse_matrix.templates.h
tests/lac/sparse_matrix_move.cc [new file with mode: 0644]
tests/lac/sparse_matrix_move.with_cxx11=on.output [new file with mode: 0644]

index 1cde275e27346db2a2794cebeefaf588e28e18e8..50ea51f48aa5ff6fe4296a58c017d98a2fe44624 100644 (file)
@@ -543,6 +543,20 @@ public:
    */
   SparseMatrix (const SparseMatrix &);
 
+#ifdef DEAL_II_WITH_CXX11
+  /**
+   * Move constructor. Construct a new sparse matrix by transferring the
+   * internal data of the matrix @p m into a new object.
+   *
+   * Move construction allows an object to be returned from a function or
+   * packed into a tuple even when the class cannot be copy-constructed.
+   *
+   * @note this constructor is only available if deal.II is configured with
+   * C++11 support.
+   */
+  SparseMatrix (SparseMatrix<number> &&m);
+#endif
+
   /**
    * Constructor. Takes the given matrix sparsity structure to represent the
    * sparsity pattern of this matrix. You can change the sparsity pattern
@@ -584,6 +598,17 @@ public:
    */
   SparseMatrix<number> &operator = (const SparseMatrix<number> &);
 
+#ifdef DEAL_II_WITH_CXX11
+  /**
+   * Move assignment operator. This operator replaces the present matrix with
+   * @p m by transferring the internal data of @p m.
+   *
+   * @note This operator is only available if deal.II is configured with C++11
+   * support.
+   */
+  SparseMatrix<number> &operator = (SparseMatrix<number> &&m);
+#endif
+
   /**
    * Copy operator: initialize the matrix with the identity matrix. This
    * operator will throw an exception if the sizes of the sparsity pattern and
index 4f145c555c72c6fe2bf52f0736e33a37d523785f..0cda654f7ad3a627659aafa26402678312249601 100644 (file)
@@ -69,6 +69,23 @@ SparseMatrix<number>::SparseMatrix (const SparseMatrix &m)
 
 
 
+#ifdef DEAL_II_WITH_CXX11
+template <typename number>
+SparseMatrix<number>::SparseMatrix (SparseMatrix<number> &&m)
+  :
+  Subscriptor (m),
+  cols(m.cols),
+  val(m.val),
+  max_len(m.max_len)
+{
+  m.cols = 0;
+  m.val = 0;
+  m.max_len = 0;
+}
+#endif
+
+
+
 template <typename number>
 SparseMatrix<number> &
 SparseMatrix<number>::operator = (const SparseMatrix<number> &m)
@@ -83,6 +100,25 @@ SparseMatrix<number>::operator = (const SparseMatrix<number> &m)
 
 
 
+#ifdef DEAL_II_WITH_CXX11
+template <typename number>
+SparseMatrix<number> &
+SparseMatrix<number>::operator = (SparseMatrix<number> &&m)
+{
+  cols = m.cols;
+  val = m.val;
+  max_len = m.max_len;
+
+  m.cols = 0;
+  m.val = 0;
+  m.max_len = 0;
+
+  return *this;
+}
+#endif
+
+
+
 template <typename number>
 SparseMatrix<number>::SparseMatrix (const SparsityPattern &c)
   :
diff --git a/tests/lac/sparse_matrix_move.cc b/tests/lac/sparse_matrix_move.cc
new file mode 100644 (file)
index 0000000..6d67eb3
--- /dev/null
@@ -0,0 +1,100 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 1998 - 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 "testmatrix.h"
+#include <deal.II/base/logstream.h>
+#include <deal.II/lac/sparse_matrix.h>
+
+//#include <fstream>
+
+void graph_laplacian(const SparsityPattern& sparsity,
+                     SparseMatrix<double>& matrix)
+{
+  matrix = 0.0;
+
+  for (SparsityPattern::const_iterator it = sparsity.begin();
+       it != sparsity.end(); ++it)
+    {
+      const auto i = (*it).row();
+      const auto j = (*it).column();
+      matrix.add(i, j, -1);
+      matrix.add(i, i, 1);
+    }
+}
+
+
+SparseMatrix<double> graph_laplacian(const SparsityPattern& sparsity)
+{
+  SparseMatrix<double> A(sparsity);
+  graph_laplacian(sparsity, A);
+
+  return A;
+}
+
+
+SparseMatrix<double>
+graph_laplacian_move_return(const SparsityPattern& sparsity)
+{
+  SparseMatrix<double> A(sparsity);
+  graph_laplacian(sparsity, A);
+
+  return std::move(A);
+}
+
+
+int main()
+{
+  initlog();
+  deallog << std::setprecision(3);
+
+  const unsigned int size = 5;
+
+  FDMatrix testproblem (size, size);
+  unsigned int dim = (size-1) * (size-1);
+
+  SparsityPattern sparsity(dim, dim, size);
+  testproblem.five_point_structure(sparsity);
+  sparsity.compress();
+
+  // Return a sparse matrix, possibly using RVO or move constructor
+  SparseMatrix<double> A = graph_laplacian(sparsity);
+  deallog << A.n_nonzero_elements() << std::endl;
+
+  Vector<double> x(dim), y(dim);
+  x = 1.0;
+  y = 1.0;
+  A.vmult(y, x);
+
+  deallog << y.l2_norm() << std::endl;
+
+  // Return a sparse matrix using the move constructor
+  SparseMatrix<double> B = graph_laplacian_move_return(sparsity);
+  deallog << B.n_nonzero_elements() << std::endl;
+  y = 1.0;
+  B.vmult(y, x);
+
+  deallog << y.l2_norm() << std::endl;
+
+  // Explicitly move a sparse matrix
+  SparseMatrix<double> C;
+  C = std::move(B);
+  deallog << C.m() << std::endl;
+  deallog << B.empty() << std::endl;
+
+  return 0;
+}
diff --git a/tests/lac/sparse_matrix_move.with_cxx11=on.output b/tests/lac/sparse_matrix_move.with_cxx11=on.output
new file mode 100644 (file)
index 0000000..7239c9c
--- /dev/null
@@ -0,0 +1,7 @@
+
+DEAL::64
+DEAL::0.00
+DEAL::64
+DEAL::0.00
+DEAL::16
+DEAL::1

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.