]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Change position of template argument in CUDAWrappers::kernel
authorDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Wed, 31 Oct 2018 18:41:09 +0000 (19:41 +0100)
committerDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Fri, 4 Jan 2019 13:26:47 +0000 (14:26 +0100)
include/deal.II/lac/cuda_kernels.h
include/deal.II/lac/cuda_kernels.templates.h
source/lac/cuda_kernels.cu

index 4059cd8a06fd79211fd3fa5a0a7e3ac63543c9b8..b375af82bde3e2c9d54c702b5f2552f72280cc07 100644 (file)
@@ -26,6 +26,8 @@
 
 #  include <deal.II/lac/cuda_atomic.h>
 
+#  include <assert.h>
+
 DEAL_II_NAMESPACE_OPEN
 
 namespace LinearAlgebra
@@ -57,9 +59,9 @@ namespace LinearAlgebra
        *
        * @ingroup CUDAWrappers
        */
+      template <typename Number>
       struct Binop_Addition
       {
-        template <typename Number>
         __device__ static inline Number
         operation(const Number a, const Number b)
         {
@@ -67,6 +69,18 @@ namespace LinearAlgebra
         }
       };
 
+      template <typename Number>
+      struct Binop_Addition<std::complex<Number>>
+      {
+        __device__ static inline std::complex<Number>
+        operation(const std::complex<Number> a, const std::complex<Number> b)
+        {
+          printf("This function is not implemented for std::complex<Number>!");
+          assert(false);
+          return {};
+        }
+      };
+
 
 
       /**
@@ -74,9 +88,9 @@ namespace LinearAlgebra
        *
        * @ingroup CUDAWrappers
        */
+      template <typename Number>
       struct Binop_Subtraction
       {
-        template <typename Number>
         __device__ static inline Number
         operation(const Number a, const Number b)
         {
@@ -91,9 +105,9 @@ namespace LinearAlgebra
        *
        * @ingroup CUDAWrappers
        */
-      template <typename Number, typename Binop>
+      template <typename Number, template <typename> class Binop>
       __global__ void
-      vector_bin_op(Number *v1, Number *v2, const size_type N);
+      vector_bin_op(Number *v1, const Number *v2, const size_type N);
 
 
 
index 3a228834535c04dc7efb0aee27b73ef788049f40..7a8185f302e22eba946dfda5e9e738102fc786d5 100644 (file)
@@ -44,9 +44,9 @@ namespace LinearAlgebra
 
 
 
-      template <typename Number, typename Binop>
+      template <typename Number, template <typename> class Binop>
       __global__ void
-      vector_bin_op(Number *v1, Number *v2, const size_type N)
+      vector_bin_op(Number *v1, const Number *v2, const size_type N)
       {
         const size_type idx_base =
           threadIdx.x + blockIdx.x * (blockDim.x * chunk_size);
@@ -54,7 +54,7 @@ namespace LinearAlgebra
           {
             const size_type idx = idx_base + i * block_size;
             if (idx < N)
-              v1[idx] = Binop::operation(v1[idx], v2[idx]);
+              v1[idx] = Binop<Number>::operation(v1[idx], v2[idx]);
           }
       }
 
index 1ec88bd11dc8ecbe4cebeb51324016cb01db63b5..49230c517cc3339e239c302e79c2be86233c1b72 100644 (file)
@@ -31,11 +31,11 @@ namespace LinearAlgebra
       vec_scale<float>(float *, const float a, const size_type);
       template __global__ void
       vector_bin_op<float, Binop_Addition>(float *         v1,
-                                           float *         v2,
+                                           const float *   v2,
                                            const size_type N);
       template __global__ void
       vector_bin_op<float, Binop_Subtraction>(float *         v1,
-                                              float *         v2,
+                                              const float *   v2,
                                               const size_type N);
       template struct ElemSum<float>;
       template struct L1Norm<float>;
@@ -126,11 +126,11 @@ namespace LinearAlgebra
       vec_scale<double>(double *, const double a, const size_type);
       template __global__ void
       vector_bin_op<double, Binop_Addition>(double *        v1,
-                                            double *        v2,
+                                            const double *  v2,
                                             const size_type N);
       template __global__ void
       vector_bin_op<double, Binop_Subtraction>(double *        v1,
-                                               double *        v2,
+                                               const double *  v2,
                                                const size_type N);
       template struct ElemSum<double>;
       template struct L1Norm<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.