From: Daniel Arndt Date: Thu, 16 Aug 2018 21:06:06 +0000 (+0200) Subject: Use unique_ptr for storing the pointer to device memory X-Git-Tag: v9.1.0-rc1~800^2~4 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=2dc10e75bde7c86cd1a8ef5ec0b60697f81a3fd5;p=dealii.git Use unique_ptr for storing the pointer to device memory --- diff --git a/include/deal.II/lac/cuda_vector.h b/include/deal.II/lac/cuda_vector.h index f976eab152..184c8c55b7 100644 --- a/include/deal.II/lac/cuda_vector.h +++ b/include/deal.II/lac/cuda_vector.h @@ -67,6 +67,11 @@ namespace LinearAlgebra */ Vector(const Vector &V); + /** + * Move constructor. + */ + Vector(Vector &&) = default; + /** * Constructor. Set dimension to @p n and initialize all elements with * zero. @@ -80,9 +85,16 @@ namespace LinearAlgebra explicit Vector(const size_type n); /** - * Destructor. + * Copy assignment operator. + */ + Vector & + operator=(const Vector &v); + + /** + * Move assignment operator. */ - ~Vector(); + Vector & + operator=(Vector &&v) = default; /** * Reinit functionality. The flag omit_zeroing_entries @@ -254,7 +266,8 @@ namespace LinearAlgebra const VectorSpaceVector &W) override; /** - * Return the pointer to the underlying array. + * Return the pointer to the underlying array. Ownership still resides + * with this class. */ Number * get_values() const; @@ -298,7 +311,7 @@ namespace LinearAlgebra /** * Pointer to the array of elements of this vector. */ - Number *val; + std::unique_ptr val; /** * Number of elements in the vector. @@ -313,7 +326,7 @@ namespace LinearAlgebra inline Number * Vector::get_values() const { - return val; + return val.get(); } diff --git a/source/lac/cuda_vector.cu b/source/lac/cuda_vector.cu index 58bdbe4972..e40b47c802 100644 --- a/source/lac/cuda_vector.cu +++ b/source/lac/cuda_vector.cu @@ -35,9 +35,35 @@ namespace LinearAlgebra using ::dealii::CUDAWrappers::block_size; using ::dealii::CUDAWrappers::chunk_size; + namespace + { + template + void + delete_device_vector(Number *device_ptr) noexcept + { + cudaError_t error_code = cudaFree(device_ptr); + (void)error_code; + AssertNothrow(error_code == cudaSuccess, + dealii::ExcCudaError(cudaGetErrorString(error_code))); + } + + template + Number * + allocate_device_vector(const std::size_t size) + { + Number * device_ptr; + cudaError_t error_code = cudaMalloc(&device_ptr, size * sizeof(Number)); + (void)error_code; + AssertCuda(error_code); + return device_ptr; + } + } // namespace + + + template Vector::Vector() - : val(nullptr) + : val(nullptr, delete_device_vector) , n_elements(0) {} @@ -45,41 +71,45 @@ namespace LinearAlgebra template Vector::Vector(const Vector &V) - : n_elements(V.n_elements) + : val(allocate_device_vector(V.n_elements), + delete_device_vector) + , n_elements(V.n_elements) { - // Allocate the memory - cudaError_t error_code = cudaMalloc(&val, n_elements * sizeof(Number)); - AssertCuda(error_code); // Copy the values. - error_code = cudaMemcpy(val, - V.val, - n_elements * sizeof(Number), - cudaMemcpyDeviceToDevice); + cudaError_t error_code = cudaMemcpy(val.get(), + V.val.get(), + n_elements * sizeof(Number), + cudaMemcpyDeviceToDevice); AssertCuda(error_code); } template - Vector::Vector(const size_type n) - : val(nullptr) - , n_elements(0) + Vector & + Vector::operator=(const Vector &V) { - reinit(n, false); + if (n_elements < V.n_elements) + reinit(V.n_elements); + + n_elements = V.n_elements; + + // Copy the values. + cudaError_t error_code = cudaMemcpy(val.get(), + V.val.get(), + n_elements * sizeof(Number), + cudaMemcpyDeviceToDevice); + AssertCuda(error_code); } template - Vector::~Vector() + Vector::Vector(const size_type n) + : val(nullptr, delete_device_vector) + , n_elements(0) { - if (val != nullptr) - { - cudaError_t error_code = cudaFree(val); - AssertCuda(error_code); - val = nullptr; - n_elements = 0; - } + reinit(n, false); } @@ -90,31 +120,15 @@ namespace LinearAlgebra { // Resize the underlying array if necessary if (n == 0) - { - if (val != nullptr) - { - cudaError_t error_code = cudaFree(val); - AssertCuda(error_code); - val = nullptr; - } - } - else - { - if ((n_elements != n) && (val != nullptr)) - { - cudaError_t error_code = cudaFree(val); - AssertCuda(error_code); - } + val.reset(); + else if (n != n_elements) + val.reset(allocate_device_vector(n)); - cudaError_t error_code = cudaMalloc(&val, n * sizeof(Number)); + // If necessary set the elements to zero + if (omit_zeroing_entries == false) + { + cudaError_t error_code = cudaMemset(val.get(), 0, n * sizeof(Number)); AssertCuda(error_code); - - // If necessary set the elements to zero - if (omit_zeroing_entries == false) - { - cudaError_t error_code = cudaMemset(val, 0, n * sizeof(Number)); - AssertCuda(error_code); - } } n_elements = n; } @@ -139,7 +153,7 @@ namespace LinearAlgebra { if (operation == VectorOperation::insert) { - cudaError_t error_code = cudaMemcpy(val, + cudaError_t error_code = cudaMemcpy(val.get(), V.begin(), n_elements * sizeof(Number), cudaMemcpyHostToDevice); @@ -164,7 +178,7 @@ namespace LinearAlgebra const int n_blocks = 1 + (n_elements - 1) / (chunk_size * block_size); kernel::vector_bin_op - <<>>(val, tmp, n_elements); + <<>>(val.get(), tmp, n_elements); // Check that the kernel was launched correctly AssertCuda(cudaGetLastError()); // Check that there was no problem during the execution of the kernel @@ -187,7 +201,8 @@ namespace LinearAlgebra Assert(s == Number(), ExcMessage("Only 0 can be assigned to a vector.")); (void)s; - cudaError_t error_code = cudaMemset(val, 0, n_elements * sizeof(Number)); + cudaError_t error_code = + cudaMemset(val.get(), 0, n_elements * sizeof(Number)); AssertCuda(error_code); return *this; @@ -202,7 +217,7 @@ namespace LinearAlgebra AssertIsFinite(factor); const int n_blocks = 1 + (n_elements - 1) / (chunk_size * block_size); kernel::vec_scale - <<>>(val, factor, n_elements); + <<>>(val.get(), factor, n_elements); // Check that the kernel was launched correctly AssertCuda(cudaGetLastError()); @@ -222,7 +237,7 @@ namespace LinearAlgebra Assert(factor != Number(0.), ExcZero()); const int n_blocks = 1 + (n_elements - 1) / (chunk_size * block_size); kernel::vec_scale - <<>>(val, 1. / factor, n_elements); + <<>>(val.get(), 1. / factor, n_elements); // Check that the kernel was launched correctly AssertCuda(cudaGetLastError()); @@ -251,7 +266,7 @@ namespace LinearAlgebra const int n_blocks = 1 + (n_elements - 1) / (chunk_size * block_size); kernel::vector_bin_op - <<>>(val, down_V.val, n_elements); + <<>>(val.get(), down_V.val.get(), n_elements); // Check that the kernel was launched correctly AssertCuda(cudaGetLastError()); @@ -280,7 +295,7 @@ namespace LinearAlgebra const int n_blocks = 1 + (n_elements - 1) / (chunk_size * block_size); kernel::vector_bin_op - <<>>(val, down_V.val, n_elements); + <<>>(val.get(), down_V.val.get(), n_elements); // Check that the kernel was launched correctly AssertCuda(cudaGetLastError()); @@ -314,8 +329,8 @@ namespace LinearAlgebra const int n_blocks = 1 + (n_elements - 1) / (chunk_size * block_size); kernel::double_vector_reduction> <<>>(result_device, - val, - down_V.val, + val.get(), + down_V.val.get(), static_cast( n_elements)); @@ -341,7 +356,8 @@ namespace LinearAlgebra { AssertIsFinite(a); const int n_blocks = 1 + (n_elements - 1) / (chunk_size * block_size); - kernel::vec_add<<>>(val, a, n_elements); + kernel::vec_add + <<>>(val.get(), a, n_elements); // Check that the kernel was launched correctly AssertCuda(cudaGetLastError()); @@ -369,7 +385,7 @@ namespace LinearAlgebra const int n_blocks = 1 + (n_elements - 1) / (chunk_size * block_size); kernel::add_aV<<>>( - val, a, down_V.val, n_elements); + val.get(), a, down_V.val.get(), n_elements); // Check that the kernel was launched correctly AssertCuda(cudaGetLastError()); @@ -411,7 +427,7 @@ namespace LinearAlgebra const int n_blocks = 1 + (n_elements - 1) / (chunk_size * block_size); kernel::add_aVbW<<>>( - val, a, down_V.val, b, down_W.val, n_elements); + val.get(), a, down_V.val.get(), b, down_W.val.get(), n_elements); // Check that the kernel was launched correctly AssertCuda(cudaGetLastError()); @@ -442,7 +458,7 @@ namespace LinearAlgebra const int n_blocks = 1 + (n_elements - 1) / (chunk_size * block_size); kernel::sadd<<>>( - s, val, a, down_V.val, n_elements); + s, val.get(), a, down_V.val.get(), n_elements); // Check that the kernel was launched correctly AssertCuda(cudaGetLastError()); @@ -468,10 +484,8 @@ namespace LinearAlgebra "Cannot scale two vectors with different numbers of elements.")); const int n_blocks = 1 + (n_elements - 1) / (chunk_size * block_size); - kernel::scale - <<>>(val, - down_scaling_factors.val, - n_elements); + kernel::scale<<>>( + val.get(), down_scaling_factors.val.get(), n_elements); // Check that the kernel was launched correctly AssertCuda(cudaGetLastError()); @@ -499,10 +513,8 @@ namespace LinearAlgebra "Cannot assign two vectors with different numbers of elements.")); const int n_blocks = 1 + (n_elements - 1) / (chunk_size * block_size); - kernel::equ<<>>(val, - a, - down_V.val, - n_elements); + kernel::equ<<>>( + val.get(), a, down_V.val.get(), n_elements); // Check that the kernel was launched correctly AssertCuda(cudaGetLastError()); @@ -533,7 +545,7 @@ namespace LinearAlgebra const int n_blocks = 1 + (n_elements - 1) / (chunk_size * block_size); kernel::reduction> <<>>(result_device, - val, + val.get(), n_elements); // Copy the result back to the host @@ -565,7 +577,7 @@ namespace LinearAlgebra const int n_blocks = 1 + (n_elements - 1) / (chunk_size * block_size); kernel::reduction> <<>>(result_device, - val, + val.get(), n_elements); // Copy the result back to the host @@ -605,7 +617,7 @@ namespace LinearAlgebra const int n_blocks = 1 + (n_elements - 1) / (chunk_size * block_size); kernel::reduction> <<>>(result_device, - val, + val.get(), n_elements); // Copy the result back to the host @@ -654,7 +666,7 @@ namespace LinearAlgebra const int n_blocks = 1 + (n_elements - 1) / (chunk_size * block_size); kernel::add_and_dot<<>>( - res_d, val, down_V.val, down_W.val, a, n_elements); + res_d, val.get(), down_V.val.get(), down_W.val.get(), a, n_elements); Number res; error_code = @@ -690,7 +702,7 @@ namespace LinearAlgebra // Copy the vector to the host std::vector cpu_val(n_elements); - Utilities::CUDA::copy_to_host(val, cpu_val); + Utilities::CUDA::copy_to_host(val.get(), cpu_val); for (unsigned int i = 0; i < n_elements; ++i) out << cpu_val[i] << std::endl; out << std::flush;