*/
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
*/
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
+#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)
+#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)
:
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// 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;
+}