// ---------------------------------------------------------------------
// $Id: paralution_sparse_matrix.h 30040 2013-07-18 17:06:48Z maier $
//
-// Copyright (C) 2013 by the deal.II authors
+// Copyright (C) 2013, 2014 by the deal.II authors
//
// This file is part of the deal.II library.
//
*/
namespace ParalutionWrappers
{
+ /**
+ * Enum on format of Paralusion Matrix: DENSE, Compressed Sparse Row,
+ * Modified Compressed Sparse Row, Block Compressed Sparse Row, COOrdinate
+ * format, DIAgonal format, ELLpack format, HYBrid format (mixed between ELL
+ * and COO).
+ */
+ enum matrix_format{DENSE, CSR, MCSR, BCSR, COO, DIA, ELL, HYB};
+
/**
* This class implements a wrapper to use the Paralution sparse matrix class
* LocalMatrix. This class is designed for use in either serial
*
* @ingroup ParalutionWrappers
* @ingroup Matrix
- * @author Bruno Turcksin, 2013
+ * @author Bruno Turcksin, 2013, 2014
*/
template <typename Number>
class SparseMatrix : public Subscriptor
SparseMatrix();
/**
- * Copy constructor. This constructor is only allowed to be called if
- * the matrix to be copied is empty. This is for the same reason as for
- * the SparsityPattern, see there fore the details.
- *
- * If you really want to copy a whole matrix, you can do so by using the
- * copy_from() function.
+ * Copy constructor. If @p copy_backend is set to false, the copied sparse
+ * matrix stays on the host/device where it is created. Otherwise the copied
+ * sparse matrix is moved to the host/device of the given vector.
*/
- //TODO
- //SparseMarix(const SparseMatrix &sparse_matrix);
+ SparseMatrix(const SparseMatrix &sparse_matrix, bool copy_backend = false);
/**
- * Constructo. Takes the given matrix sparsity structure to represent
+ * Constructor. Takes the given matrix sparsity structure to represent
* the sparsity pattern of this matrix. You can change the sparsity
* pattern later on by calling the reinit(const SparsityPattern&)
* function.
/**
* This function converts the underlying SparseMatrix to
- * Paralution::LocalMatrix. This function frees the SparseMatrix.
+ * Paralution::LocalMatrix. This function frees the dealii::SparseMatrix.
*/
void convert_to_paralution_csr();
//@}
* of dimension $m \times n$.
*/
size_type n() const;
+
+ /**
+ * This function returns true if the underlying matrix is a
+ * paralution::LocalMatrix.
+ */
+ bool is_paralution_matrix() const;
+
+ /**
+ * Return the format of the Paralution::LocalMatrix.
+ */
+ matrix_format get_matrix_format() const;
//@}
/**
* @name 3: Modifying entries
/**
* Same function as before, but now including the possibility to use
- * rectaangular full_matrices and different local-to-global indexing on
+ * rectangular full_matrices and different local-to-global indexing on
* rows and columns, respectively.
*
* This function only works on the host.
const bool elide_zero_values = true,
const bool col_indices_are_sorted = false);
+ /**
+ * Transpose the matrix. This function can only be called if the underlying
+ * is a Paralution::LocalMatrix.
+ */
+ void transpose();
+
+ /**
+ * Copy @p sparse_matrix.
+ */
+ void copy_from(const SparseMatrix &sparse_matrix);
+
+ /**
+ * Copy @p sparse_matrix. If the underlying matrix is a
+ * Paralution::LocalMatrix, the functions returns immediately and performs
+ * the asynchronous transder in the background.
+ */
+ void copy_from_async(const SparseMatrix &sparse_matrix);
+
//@}
/**
- * @name 4: Access to underlying Paralution data and move data to
+ * @name 4: Access to underlying data and move data to
* accelerator/host
*/
/**
void move_to_host();
/**
- * Return a constant reference to the underlying Paralution LocalMatrix
+ * Move the SparseMatrix to the accelerator. The function returns
+ * immediately and performs the asynchronous transfer in the background.
+ * This function should only be called after convert_to_paralution_csr.
+ */
+ void move_to_accelerator_async();
+
+ /**
+ * Move the SparseMatrix to the host. The function returns immediately and
+ * performs the asynchronous transfer in the background. This function
+ * should only be called after convert_to_paralution_csr.
+ */
+ void move_to_host_async();
+
+ /**
+ * Synchronize the code when move_to_host_async or move_to_accelerator_async
+ * is used.
+ */
+ void sync();
+
+ /**
+ * Return a constant reference to the underlying dealii::SparseMatrix.
+ */
+ ::dealii::SparseMatrix<Number> const &dealii_matrix() const;
+
+ /**
+ * Return a reference to the underlying dealii::SparseMatrix.
+ */
+ ::dealii::SparseMatrix<Number>& dealii_matrix();
+
+ /**
+ * Return a constant reference to the underlying Paralution::LocalMatrix
* data.
*/
paralution::LocalMatrix<Number> const ¶lution_matrix() const;
/**
- * Return a reference to the underlying Paralution LocalMatrix data.
+ * Return a reference to the underlying Paralution::LocalMatrix data.
*/
paralution::LocalMatrix<Number>& paralution_matrix();
//@}
private :
/**
- * This flag is true if local_matrix is used. It becomes true after
- * convert_to_paralution_csr has been called.
+ * This flag is true if @p local_matrix is used. It becomes true after
+ * convert_to_paralution_csr has been called, i.e. after the @p
+ * sparse_matrix as be freed and the paralution::LocalMatrix has been
+ * created.
*/
bool is_local_matrix;
+ template <typename Number>
+ inline bool SparseMatrix<Number>::is_paralution_matrix() const
+ {
+ return is_local_matrix;
+ }
+
+
+
+ template <typename Number>
+ inline matrix_format SparseMatrix<Number>::get_matrix_format() const
+ {
+ return local_matrix.get_format();
+ }
+
+
+
template <typename Number>
inline void SparseMatrix<Number>::set(const size_type i,
const size_type j,
+ template <typename Number>
+ inline void SparseMatrix<Number>::move_to_accelerator_async()
+ {
+ local_matrix.MoveToAcceleratorAsync();
+ }
+
+
+
+ template <typename Number>
+ inline void SparseMatrix<Number>::move_to_host_async()
+ {
+ local_matrix.MoveToHostAsync();
+ }
+
+
+
+ template <typename Number>
+ inline void SparseMatrix<Number>::sync()
+ {
+ local_matrix.Sync();
+ }
+
+
+
+ template <typename Number>
+ inline ::dealii::SparseMatrix<Number> const &SparseMatrix<Number>::dealii_matrix() const
+ {
+ return sparse_matrix;
+ }
+
+
+
+ template <typename Number>
+ inline ::dealii::SparseMatrix<Number>& SparseMatrix<Number>::dealii_matrix()
+ {
+ return sparse_matrix;
+ }
+
+
+
template <typename Number>
inline paralution::LocalMatrix<Number> const &SparseMatrix<Number>::paralution_matrix() const
{
// ---------------------------------------------------------------------
// $Id: paralution_sparse_matrix.cc 31567 2013-11-06 18:01:36Z turcksin $
//
-// Copyright (C) 2013 by the deal.II authors
+// Copyright (C) 2013, 2014 by the deal.II authors
//
// This file is part of the deal.II library.
//
namespace ParalutionWrappers
{
+ template <typename Number>
+ SparseMatrix<Number>::SparseMatrix(const SparseMatrix &sparse_matrix, bool copy_backend)
+ {
+ is_local_matrix = sparse_matrix.is_paralution_matrix();
+ // This is costly but copy_from is used to have a consistent behavior between
+ // sparse_matrix and local_matrix.
+ if (is_local_matrix==false)
+ sparse_matrix.copy_from(sparse_matrix.sparse_matrix());
+ else
+ {
+ if (copy_backend==false)
+ local_matrix.CopyFrom(sparse_matrix.paralution_matrix());
+ else
+ local_matrix.CloneFrom(sparse_matrix.paralution_matrix());
+ }
+ }
+
+
+
template <typename Number>
void SparseMatrix<Number>::convert_to_paralution_csr()
{
sparse_matrix.clear();
is_local_matrix = true;
}
+
+
+
+ template <typename Number>
+ void SparseMatrix<Number>::transpose()
+ {
+ if (is_local_matrix==false)
+ {
+ AssertThrow(false,
+ ExcMessage("Transpose cannot be used if the underlying matrix is a dealii::SparseMatrix."));
+ }
+ else
+ {
+ local_matrix.Transpose();
+ }
+ }
+
+
+
+ template <typename Number>
+ void SparseMatrix<Number>::copy_from(const SparseMatrix &sparse_matrix)
+ {
+ is_local_matrix = sparse_matrix.is_paralution_matrix();
+ if (is_local_matrix==false)
+ sparse_matrix.copy_from(sparse_matrix.sparse_matrix());
+ else
+ local_matrix.CopyFrom(sparse_matrix.paralution_matrix());
+ }
+
+
+
+ template <typename Number>
+ void SparseMatrix<Number>::copy_from_async(const SparseMatrix &sparse_matrix)
+ {
+ is_local_matrix = sparse_matrix.is_paralution_matrix();
+ if (is_local_matrix==false)
+ sparse_matrix.copy_from(sparse_matrix.sparse_matrix());
+ else
+ local_matrix.CopyFromAsync(sparse_matrix.paralution_matrix());
+ }
}
// Explicit instantiations