]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add new cuda kernel function to execute a binary operation on a subset of the
authorBruno Turcksin <bruno.turcksin@gmail.com>
Fri, 8 Feb 2019 14:47:25 +0000 (14:47 +0000)
committerBruno Turcksin <bruno.turcksin@gmail.com>
Tue, 12 Feb 2019 01:32:31 +0000 (01:32 +0000)
vector

include/deal.II/lac/cuda_kernels.h
include/deal.II/lac/cuda_kernels.templates.h
source/lac/cuda_kernels.cu

index 963a9461702a8b929fb4b9cb1ef35f207df738ef..c3470a5749d983977380c6faa676e14e131094ea 100644 (file)
@@ -181,6 +181,23 @@ namespace LinearAlgebra
 
 
 
+      /**
+       * Apply the functor @tparam Binop to the elements of @p v1 that have
+       * indices in @p mask and @p v2. The size of @p mask should be greater
+       * than the size of @p v1. @p mask and @p v2 should have the same size @p
+       * N.
+       *
+       * @ingroup CUDAWrappers
+       */
+      template <typename Number, template <typename> class Binop>
+      __global__ void
+      masked_vector_bin_op(const unsigned int *mask,
+                           Number *            v1,
+                           const Number *      v2,
+                           const size_type     N);
+
+
+
       /**
        * Structure implementing the functions used to add elements when
        * using a reduction.
index cd522131ba9c1c4bb822b0511e3196a09dc6c11f..cadeb4d11a39eda8708dc1097569def637bb4f27 100644 (file)
@@ -60,6 +60,24 @@ namespace LinearAlgebra
 
 
 
+      template <typename Number, template <typename> class Binop>
+      __global__ void
+      masked_vector_bin_op(const unsigned int *mask,
+                           Number *            v1,
+                           const Number *      v2,
+                           const size_type     N)
+      {
+        const size_type idx_base =
+          threadIdx.x + blockIdx.x * (blockDim.x * chunk_size);
+        for (unsigned int i = 0; i < chunk_size; ++i)
+          {
+            const size_type idx = idx_base + i * block_size;
+            if (idx < N)
+              v1[mask[idx]] = Binop<Number>::operation(v1[mask[idx]], v2[idx]);
+          }
+      }
+
+
       template <typename Number>
       __device__ Number
                  ElemSum<Number>::reduction_op(const Number a, const Number b)
index a1c83d59352d4de6c4cd8d799072b84b8932334f..c7bfb9b020e3d3ace05c629bdb67488f80a8fd26 100644 (file)
@@ -37,6 +37,16 @@ namespace LinearAlgebra
       vector_bin_op<float, Binop_Subtraction>(float *         v1,
                                               const float *   v2,
                                               const size_type N);
+      template __global__ void
+      masked_vector_bin_op<float, Binop_Addition>(const unsigned int *mask,
+                                                  float *             v1,
+                                                  const float *       v2,
+                                                  const size_type     N);
+      template __global__ void
+      masked_vector_bin_op<float, Binop_Subtraction>(const unsigned int *mask,
+                                                     float *             v1,
+                                                     const float *       v2,
+                                                     const size_type     N);
       template struct ElemSum<float>;
       template struct L1Norm<float>;
       template struct LInfty<float>;
@@ -137,6 +147,16 @@ namespace LinearAlgebra
       vector_bin_op<double, Binop_Subtraction>(double *        v1,
                                                const double *  v2,
                                                const size_type N);
+      template __global__ void
+      masked_vector_bin_op<double, Binop_Addition>(const unsigned int *mask,
+                                                   double *            v1,
+                                                   const double *      v2,
+                                                   const size_type     N);
+      template __global__ void
+      masked_vector_bin_op<double, Binop_Subtraction>(const unsigned int *mask,
+                                                      double *            v1,
+                                                      const double *      v2,
+                                                      const size_type     N);
       template struct ElemSum<double>;
       template struct L1Norm<double>;
       template struct LInfty<double>;

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.