+ /**
+ * 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.
+ 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)
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>;
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>;