From 3500065bb028ae9fc388ba84817d29537f3c3046 Mon Sep 17 00:00:00 2001 From: Bruno Turcksin Date: Sun, 17 Dec 2017 17:27:42 -0500 Subject: [PATCH] Add move constructor for TrilinosWrappers::SparseMatrix --- include/deal.II/lac/trilinos_sparse_matrix.h | 25 +++++--- source/lac/trilinos_sparse_matrix.cc | 15 +++++ tests/trilinos/sparse_matrix_08.cc | 60 ++++++++++++++++++++ tests/trilinos/sparse_matrix_08.output | 2 + 4 files changed, 93 insertions(+), 9 deletions(-) create mode 100644 tests/trilinos/sparse_matrix_08.cc create mode 100644 tests/trilinos/sparse_matrix_08.output diff --git a/include/deal.II/lac/trilinos_sparse_matrix.h b/include/deal.II/lac/trilinos_sparse_matrix.h index c3e3676526..0eac855c31 100644 --- a/include/deal.II/lac/trilinos_sparse_matrix.h +++ b/include/deal.II/lac/trilinos_sparse_matrix.h @@ -567,6 +567,22 @@ namespace TrilinosWrappers */ SparseMatrix (const SparsityPattern &InputSparsityPattern); + /** + * Move constructor. Create a new sparse matrix by stealing the internal + * data. + */ + SparseMatrix (SparseMatrix &&other); + + /** + * Copy constructor is deleted. + */ + SparseMatrix (const SparseMatrix &) = delete; + + /** + * operator= is deleted. + */ + SparseMatrix &operator = (const SparseMatrix &) = delete; + /** * Destructor. Made virtual so that one can use pointers to this class. */ @@ -1990,15 +2006,6 @@ namespace TrilinosWrappers private: - /** - * Copy constructor is disabled. - */ - SparseMatrix (const SparseMatrix &); - /** - * operator= is disabled. - */ - SparseMatrix &operator = (const SparseMatrix &); - /** * Pointer to the user-supplied Epetra Trilinos mapping of the matrix * columns that assigns parts of the matrix to the individual processes. diff --git a/source/lac/trilinos_sparse_matrix.cc b/source/lac/trilinos_sparse_matrix.cc index 2daa3598b2..8e502c9966 100644 --- a/source/lac/trilinos_sparse_matrix.cc +++ b/source/lac/trilinos_sparse_matrix.cc @@ -356,6 +356,21 @@ namespace TrilinosWrappers + SparseMatrix::SparseMatrix (SparseMatrix &&other) + : + column_space_map(std::move(other.column_space_map)), + matrix(std::move(other.matrix)), + nonlocal_matrix(std::move(other.nonlocal_matrix)), + nonlocal_matrix_exporter(std::move(other.nonlocal_matrix_exporter)), + last_action(other.last_action), + compressed(other.compressed) + { + other.last_action = Zero; + other.compressed = false; + } + + + void SparseMatrix::copy_from (const SparseMatrix &rhs) { diff --git a/tests/trilinos/sparse_matrix_08.cc b/tests/trilinos/sparse_matrix_08.cc new file mode 100644 index 0000000000..ab1ee3b1fe --- /dev/null +++ b/tests/trilinos/sparse_matrix_08.cc @@ -0,0 +1,60 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2017 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 move constructor + +#include "../tests.h" +#include +#include +#include +#include + +int main (int argc,char **argv) +{ + initlog(); + + Utilities::MPI::MPI_InitFinalize mpi_initialization (argc, argv, testing_max_num_threads()); + + IndexSet partitioning (3); + + partitioning.add_range(0, 3); + + // Add element (2,1) to the matrix + TrilinosWrappers::SparsityPattern sp (partitioning); + sp.add (2,1); + sp.compress(); + + TrilinosWrappers::SparseMatrix A (sp); + A.add (2, 1, 2.0); + A.compress(VectorOperation::add); + + TrilinosWrappers::SparseMatrix B(std::move(A)); + + for (unsigned int i=0; i<3; ++i) + for (unsigned int j=0; j<3; ++j) + { + if ((i == 2) && (j == 1)) + { + AssertThrow(B.el(i,j) == 2, ExcInternalError()); + } + else + { + AssertThrow(B.el(i,j) == 0, ExcInternalError()); + } + } + + deallog << "OK" << std::endl; +} diff --git a/tests/trilinos/sparse_matrix_08.output b/tests/trilinos/sparse_matrix_08.output new file mode 100644 index 0000000000..0fd8fc12f0 --- /dev/null +++ b/tests/trilinos/sparse_matrix_08.output @@ -0,0 +1,2 @@ + +DEAL::OK -- 2.39.5