--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2018 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.
+//
+// ---------------------------------------------------------------------
+
+#ifndef dealii_cuda_h
+#define dealii_cuda_h
+
+#include <deal.II/base/config.h>
+
+#ifdef DEAL_II_WITH_CUDA
+
+#include <cusparse.h>
+
+DEAL_II_NAMESPACE_OPEN
+namespace Utilities
+{
+ /**
+ * A namespace for utility structures for CUDA.
+ */
+ namespace CUDA
+ {
+ /**
+ * This structure creates, stores, and destroys the handles of the different
+ * CUDA libraries used inside deal.II.
+ */
+ struct Handle
+ {
+ /**
+ * Constructor. Create the handles for the different libraries.
+ */
+ Handle();
+
+ /**
+ * Copy constructor is deleted.
+ */
+ Handle(Handle const &) = delete;
+
+ /**
+ * Destructor. Destroy the handles and free all the memory allocated by
+ * GrowingVectorMemory.
+ */
+ ~Handle();
+
+ /**
+ * Handle to the cuSPARSE library.
+ */
+ cusparseHandle_t cusparse_handle;
+ };
+ }
+}
+
+DEAL_II_NAMESPACE_CLOSE
+
+#endif
+
+#endif
#include <deal.II/base/subscriptor.h>
#ifdef DEAL_II_WITH_CUDA
+#include <deal.II/base/cuda.h>
#include <deal.II/lac/sparse_matrix.h>
#include <deal.II/lac/cuda_vector.h>
#include <cusparse.h>
SparseMatrix();
/**
- * Constructor. Takes a cuSPARSE handle and a sparse matrix on the host.
+ * Constructor. Takes a Utilities::CUDA::Handle and a sparse matrix on the host.
* The sparse matrix on the host is copied on the device and the elements
* are reordered according to the format supported by cuSPARSE.
*/
- SparseMatrix(cusparseHandle_t handle,
+ SparseMatrix(Utilities::CUDA::Handle &handle,
const ::dealii::SparseMatrix<Number> &sparse_matrix_host);
/**
* to the device and the elementes are reordered according to the format
* supported by cuSPARSE.
*/
- void reinit(cusparseHandle_t handle,
+ void reinit(Utilities::CUDA::Handle &handle,
const ::dealii::SparseMatrix<Number> &sparse_matrix_host);
//@}
private:
/**
- * cuSPARSE used to call cuSPARSE function. The cuSPARSE handle needs to
- * be mutable to be called in a const function.
+ * cuSPARSE handle used to call cuSPARSE functions. The cuSPARSE handle needs
+ * to be mutable to be called in a const function.
*/
mutable cusparseHandle_t cusparse_handle;
symmetric_tensor.cc
)
+# Add CUDA wrapper files
+IF(DEAL_II_WITH_CUDA)
+ SET(_separate_src
+ ${_separate_src}
+ cuda.cu
+ )
+ENDIF()
+
# determined by profiling
SET(_n_includes_per_unity_file 29)
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2018 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 <deal.II/base/cuda.h>
+
+#include <deal.II/base/exceptions.h>
+#include <deal.II/lac/cuda_vector.h>
+#include <deal.II/lac/vector_memory.h>
+
+DEAL_II_NAMESPACE_OPEN
+
+namespace Utilities
+{
+ namespace CUDA
+ {
+ Handle::Handle()
+ {
+ cusparseStatus_t cusparse_error_code = cusparseCreate(&cusparse_handle);
+ AssertCusparse(cusparse_error_code);
+ }
+
+
+
+ Handle::~Handle()
+ {
+ dealii::GrowingVectorMemory<LinearAlgebra::CUDAWrappers::Vector<float>>
+ ::release_unused_memory();
+ dealii::GrowingVectorMemory<LinearAlgebra::CUDAWrappers::Vector<double>>
+ ::release_unused_memory();
+
+ cusparseStatus_t cusparse_error_code = cusparseDestroy(cusparse_handle);
+ AssertCusparse(cusparse_error_code);
+ }
+ }
+}
+
+DEAL_II_NAMESPACE_CLOSE
template <typename Number>
- SparseMatrix<Number>::SparseMatrix(cusparseHandle_t handle,
+ SparseMatrix<Number>::SparseMatrix(Utilities::CUDA::Handle &handle,
const ::dealii::SparseMatrix<Number> &sparse_matrix_host)
:
val_dev(nullptr),
template <typename Number>
- void SparseMatrix<Number>::reinit(cusparseHandle_t handle,
+ void SparseMatrix<Number>::reinit(Utilities::CUDA::Handle &handle,
const ::dealii::SparseMatrix<Number> &sparse_matrix_host)
{
- cusparse_handle = handle;
+ cusparse_handle = handle.cusparse_handle;
nnz = sparse_matrix_host.n_nonzero_elements();
n_rows = sparse_matrix_host.m();
n_cols = sparse_matrix_host.n();
#include "../tests.h"
#include "../testmatrix.h"
+#include <deal.II/base/cuda.h>
#include <deal.II/base/exceptions.h>
#include <deal.II/lac/cuda_sparse_matrix.h>
#include <deal.II/lac/precondition.h>
#include <deal.II/lac/solver_control.h>
#include <deal.II/lac/vector.h>
-void test(cusparseHandle_t cusparse_handle)
+void test(Utilities::CUDA::Handle &cuda_handle)
{
// Build the sparse matrix on the host
const unsigned int problem_size = 10;
cg_host.solve(A, sol_host, rhs_host, prec_no);
// Solve on the device
- CUDAWrappers::SparseMatrix<double> A_dev(cusparse_handle, A);
+ CUDAWrappers::SparseMatrix<double> A_dev(cuda_handle, A);
LinearAlgebra::CUDAWrappers::Vector<double> sol_dev(size);
LinearAlgebra::CUDAWrappers::Vector<double> rhs_dev(size);
LinearAlgebra::ReadWriteVector<double> rw_vector(size);
initlog();
deallog.depth_console(0);
- cusparseHandle_t cusparse_handle;
- cusparseStatus_t cusparse_error_code = cusparseCreate(&cusparse_handle);
- AssertCusparse(cusparse_error_code);
-
- test(cusparse_handle);
+ Utilities::CUDA::Handle cuda_handle;
+ test(cuda_handle);
GrowingVectorMemory<LinearAlgebra::CUDAWrappers::Vector<double>>::release_unused_memory();
- cusparse_error_code = cusparseDestroy(cusparse_handle);
- AssertCusparse(cusparse_error_code);
-
deallog << "OK" <<std::endl;
return 0;
#include "../tests.h"
#include "../testmatrix.h"
+#include <deal.II/base/cuda.h>
#include <deal.II/lac/cuda_sparse_matrix.h>
#include <deal.II/lac/read_write_vector.h>
#include <deal.II/lac/vector.h>
AssertThrow(std::abs(a[i] - b[i]) < 1e-15, ExcInternalError());
}
-void test(cusparseHandle_t cusparse_handle)
+void test(Utilities::CUDA::Handle &cuda_handle)
{
// Build the sparse matrix on the host
const unsigned int size = 10;
testproblem.five_point(A, true);
// Create the sparse matrix on the device
- CUDAWrappers::SparseMatrix<double> A_dev(cusparse_handle, A);
+ CUDAWrappers::SparseMatrix<double> A_dev(cuda_handle, A);
check_matrix(A, A_dev);
AssertThrow(A.m() == A_dev.m(), ExcInternalError());
if (i<vector_size-1)
B.set(i, i+1, 1);
}
- CUDAWrappers::SparseMatrix<double> B_dev(cusparse_handle, B);
+ CUDAWrappers::SparseMatrix<double> B_dev(cuda_handle, B);
value = B.l1_norm();
value_host = B_dev.l1_norm();
AssertThrow(std::abs(value-value_host) < 1e-15, ExcInternalError());
initlog();
deallog.depth_console(0);
- cusparseHandle_t cusparse_handle;
- cusparseStatus_t cusparse_error_code = cusparseCreate(&cusparse_handle);
- AssertCusparse(cusparse_error_code);
+ Utilities::CUDA::Handle cuda_handle;
- test(cusparse_handle);
-
- cusparse_error_code = cusparseDestroy(cusparse_handle);
- AssertCusparse(cusparse_error_code);
+ test(cuda_handle);
deallog << "OK" <<std::endl;