]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Encapsulate all the cuda handles in a struct 6269/head
authorBruno Turcksin <bruno.turcksin@gmail.com>
Fri, 13 Apr 2018 20:44:58 +0000 (16:44 -0400)
committerBruno Turcksin <bruno.turcksin@gmail.com>
Wed, 18 Apr 2018 21:04:56 +0000 (21:04 +0000)
include/deal.II/base/cuda.h [new file with mode: 0644]
include/deal.II/lac/cuda_sparse_matrix.h
source/base/CMakeLists.txt
source/base/cuda.cu [new file with mode: 0644]
source/lac/cuda_sparse_matrix.cu
tests/cuda/solver_01.cu
tests/cuda/sparse_matrix_01.cu

diff --git a/include/deal.II/base/cuda.h b/include/deal.II/base/cuda.h
new file mode 100644 (file)
index 0000000..c2594db
--- /dev/null
@@ -0,0 +1,67 @@
+// ---------------------------------------------------------------------
+//
+// 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
index 721363c82f094096e60d33c5fea05467e0106fdb..2f1a6a99fdeaf0b77f2d9aeb8c4d08e5b9a2921f 100644 (file)
@@ -20,6 +20,7 @@
 #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>
@@ -74,11 +75,11 @@ namespace CUDAWrappers
     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);
 
     /**
@@ -102,7 +103,7 @@ namespace CUDAWrappers
      * 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);
     //@}
 
@@ -251,8 +252,8 @@ namespace CUDAWrappers
 
   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;
 
index 5e80792a1a60580afa07d97d06ccf7fec8c756a8..a95b341c0d89f7e2a7cd8ab2ea2fa7f52dd7d6e7 100644 (file)
@@ -88,6 +88,14 @@ SET(_separate_src
   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)
 
diff --git a/source/base/cuda.cu b/source/base/cuda.cu
new file mode 100644 (file)
index 0000000..d89c68e
--- /dev/null
@@ -0,0 +1,49 @@
+// ---------------------------------------------------------------------
+//
+// 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
index 467d5a51f245998778edddfb4dd8c57d6eb02bd3..453330f76ea158baf4353528536aa8bc1a2fde6e 100644 (file)
@@ -134,7 +134,7 @@ namespace CUDAWrappers
 
 
   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),
@@ -208,10 +208,10 @@ namespace CUDAWrappers
 
 
   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();
index 48940f0a626e5af71ce7e9abbd638367827003e6..647aa438ec80e5f4b04bdc197d8b7f660cc8aa77 100644 (file)
@@ -18,6 +18,7 @@
 #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>
@@ -26,7 +27,7 @@
 #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;
@@ -50,7 +51,7 @@ void test(cusparseHandle_t cusparse_handle)
   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);
@@ -71,17 +72,11 @@ int main()
   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;
index e1bc483cc743d77eb1012b956d056ef5cc102ea4..20178acb45e68b5dd7f8ba8e0e571c5fbd487cd2 100644 (file)
@@ -18,6 +18,7 @@
 #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>
@@ -65,7 +66,7 @@ void check_vector(Vector<double> const &a,
     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;
@@ -79,7 +80,7 @@ void test(cusparseHandle_t cusparse_handle)
   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());
@@ -190,7 +191,7 @@ void test(cusparseHandle_t cusparse_handle)
       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());
@@ -201,14 +202,9 @@ int main()
   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;
 

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.