Loading...
Searching...
No Matches
TensorReductionCuda.h
1// This file is part of Eigen, a lightweight C++ template library
2// for linear algebra.
3//
4// Copyright (C) 2014 Benoit Steiner <benoit.steiner.goog@gmail.com>
5//
6// This Source Code Form is subject to the terms of the Mozilla
7// Public License v. 2.0. If a copy of the MPL was not distributed
8// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9
10#ifndef EIGEN_CXX11_TENSOR_TENSOR_REDUCTION_CUDA_H
11#define EIGEN_CXX11_TENSOR_TENSOR_REDUCTION_CUDA_H
12
13namespace Eigen {
14namespace internal {
15
16
17#if defined(EIGEN_USE_GPU) && defined(__CUDACC__)
18// Full reducers for GPU, don't vectorize for now
19
20// Reducer function that enables multiple cuda thread to safely accumulate at the same
21// output address. It basically reads the current value of the output variable, and
22// attempts to update it with the new value. If in the meantime another cuda thread
23// updated the content of the output address it will try again.
24template <typename T, typename R>
25__device__ EIGEN_ALWAYS_INLINE void atomicReduce(T* output, T accum, R& reducer) {
26#if __CUDA_ARCH__ >= 300
27 if (sizeof(T) == 4)
28 {
29 unsigned int oldval = *reinterpret_cast<unsigned int*>(output);
30 unsigned int newval = oldval;
31 reducer.reduce(accum, reinterpret_cast<T*>(&newval));
32 if (newval == oldval) {
33 return;
34 }
35 unsigned int readback;
36 while ((readback = atomicCAS((unsigned int*)output, oldval, newval)) != oldval) {
37 oldval = readback;
38 newval = oldval;
39 reducer.reduce(accum, reinterpret_cast<T*>(&newval));
40 if (newval == oldval) {
41 return;
42 }
43 }
44 }
45 else if (sizeof(T) == 8) {
46 unsigned long long oldval = *reinterpret_cast<unsigned long long*>(output);
47 unsigned long long newval = oldval;
48 reducer.reduce(accum, reinterpret_cast<T*>(&newval));
49 if (newval == oldval) {
50 return;
51 }
52 unsigned long long readback;
53 while ((readback = atomicCAS((unsigned long long*)output, oldval, newval)) != oldval) {
54 oldval = readback;
55 newval = oldval;
56 reducer.reduce(accum, reinterpret_cast<T*>(&newval));
57 if (newval == oldval) {
58 return;
59 }
60 }
61 }
62 else {
63 assert(0 && "Wordsize not supported");
64 }
65#else
66 EIGEN_UNUSED_VARIABLE(output);
67 EIGEN_UNUSED_VARIABLE(accum);
68 EIGEN_UNUSED_VARIABLE(reducer);
69 assert(0 && "Shouldn't be called on unsupported device");
70#endif
71}
72
73// We extend atomicExch to support extra data types
74template <typename Type>
75__device__ inline Type atomicExchCustom(Type* address, Type val) {
76 return atomicExch(address, val);
77}
78
79template <>
80__device__ inline double atomicExchCustom(double* address, double val) {
81 unsigned long long int* address_as_ull = reinterpret_cast<unsigned long long int*>(address);
82 return __longlong_as_double(atomicExch(address_as_ull, __double_as_longlong(val)));
83}
84
85#ifdef EIGEN_HAS_CUDA_FP16
86template <template <typename T> class R>
87__device__ inline void atomicReduce(half2* output, half2 accum, R<half>& reducer) {
88 unsigned int oldval = *reinterpret_cast<unsigned int*>(output);
89 unsigned int newval = oldval;
90 reducer.reducePacket(accum, reinterpret_cast<half2*>(&newval));
91 if (newval == oldval) {
92 return;
93 }
94 unsigned int readback;
95 while ((readback = atomicCAS((unsigned int*)output, oldval, newval)) != oldval) {
96 oldval = readback;
97 newval = oldval;
98 reducer.reducePacket(accum, reinterpret_cast<half2*>(&newval));
99 if (newval == oldval) {
100 return;
101 }
102 }
103}
104#endif
105
106template <>
107__device__ inline void atomicReduce(float* output, float accum, SumReducer<float>&) {
108#if __CUDA_ARCH__ >= 300
109 atomicAdd(output, accum);
110#else
111 EIGEN_UNUSED_VARIABLE(output);
112 EIGEN_UNUSED_VARIABLE(accum);
113 assert(0 && "Shouldn't be called on unsupported device");
114#endif
115}
116
117
118template <typename CoeffType, typename Index>
119__global__ void ReductionInitKernel(const CoeffType val, Index num_preserved_coeffs, CoeffType* output) {
120 const Index thread_id = blockIdx.x * blockDim.x + threadIdx.x;
121 const Index num_threads = blockDim.x * gridDim.x;
122 for (Index i = thread_id; i < num_preserved_coeffs; i += num_threads) {
123 output[i] = val;
124 }
125}
126
127
128template <int BlockSize, int NumPerThread, typename Self,
129 typename Reducer, typename Index>
130__global__ void FullReductionKernel(Reducer reducer, const Self input, Index num_coeffs,
131 typename Self::CoeffReturnType* output, unsigned int* semaphore) {
132#if __CUDA_ARCH__ >= 300
133 // Initialize the output value
134 const Index first_index = blockIdx.x * BlockSize * NumPerThread + threadIdx.x;
135 if (gridDim.x == 1) {
136 if (first_index == 0) {
137 *output = reducer.initialize();
138 }
139 }
140 else {
141 if (threadIdx.x == 0) {
142 unsigned int block = atomicCAS(semaphore, 0u, 1u);
143 if (block == 0) {
144 // We're the first block to run, initialize the output value
145 atomicExchCustom(output, reducer.initialize());
146 __threadfence();
147 atomicExch(semaphore, 2u);
148 }
149 else {
150 // Wait for the first block to initialize the output value.
151 // Use atomicCAS here to ensure that the reads aren't cached
152 unsigned int val;
153 do {
154 val = atomicCAS(semaphore, 2u, 2u);
155 }
156 while (val < 2u);
157 }
158 }
159 }
160
161 __syncthreads();
162
163 eigen_assert(gridDim.x == 1 || *semaphore >= 2u);
164
165 typename Self::CoeffReturnType accum = reducer.initialize();
166 Index max_iter = numext::mini<Index>(num_coeffs - first_index, NumPerThread*BlockSize);
167 for (Index i = 0; i < max_iter; i+=BlockSize) {
168 const Index index = first_index + i;
169 eigen_assert(index < num_coeffs);
170 typename Self::CoeffReturnType val = input.m_impl.coeff(index);
171 reducer.reduce(val, &accum);
172 }
173
174#pragma unroll
175 for (int offset = warpSize/2; offset > 0; offset /= 2) {
176#if defined(EIGEN_CUDA_SDK_VER) && EIGEN_CUDA_SDK_VER < 90000
177
178 reducer.reduce(__shfl_down(accum, offset, warpSize), &accum);
179 #else
180 reducer.reduce(__shfl_down_sync(0xFFFFFFFF, accum, offset, warpSize), &accum);
181#endif
182 }
183
184 if ((threadIdx.x & (warpSize - 1)) == 0) {
185 atomicReduce(output, accum, reducer);
186 }
187
188 if (gridDim.x > 1 && threadIdx.x == 0) {
189 // Let the last block reset the semaphore
190 atomicInc(semaphore, gridDim.x + 1);
191 }
192#else
193 EIGEN_UNUSED_VARIABLE(reducer);
194 EIGEN_UNUSED_VARIABLE(input);
195 EIGEN_UNUSED_VARIABLE(num_coeffs);
196 EIGEN_UNUSED_VARIABLE(output);
197 EIGEN_UNUSED_VARIABLE(semaphore);
198 assert(0 && "Shouldn't be called on unsupported device");
199#endif
200}
201
202
203#ifdef EIGEN_HAS_CUDA_FP16
204template <typename Self,
205 typename Reducer, typename Index>
206__global__ void ReductionInitFullReduxKernelHalfFloat(Reducer reducer, const Self input, Index num_coeffs, half2* scratch) {
207#if defined(__CUDA_ARCH__) && __CUDA_ARCH__ >= 300
208 eigen_assert(blockDim.x == 1);
209 eigen_assert(gridDim.x == 1);
210 typedef packet_traits<Eigen::half>::type packet_type;
211 Index packet_remainder = num_coeffs % Index(unpacket_traits<packet_type>::size);
212 if (packet_remainder != 0) {
213 half2* h2scratch = reinterpret_cast<half2*>(scratch);
214 for (Index i = num_coeffs - packet_remainder; i + 2 <= num_coeffs; i += 2) {
215 *h2scratch = __halves2half2(input.coeff(i), input.coeff(i + 1));
216 h2scratch++;
217 }
218 if ((num_coeffs & 1) != 0) {
219 half lastCoeff = input.coeff(num_coeffs - 1);
220 *h2scratch = __halves2half2(lastCoeff, reducer.initialize());
221 }
222 } else {
223 packet_type reduce = reducer.template initializePacket<packet_type>();
224 internal::pstoreu(scratch, reduce);
225 }
226#else
227 EIGEN_UNUSED_VARIABLE(input);
228 EIGEN_UNUSED_VARIABLE(reducer);
229 EIGEN_UNUSED_VARIABLE(num_coeffs);
230 EIGEN_UNUSED_VARIABLE(scratch);
231#endif
232}
233
234template <typename Self,
235 typename Reducer, typename Index>
236__global__ void ReductionInitKernelHalfFloat(Reducer reducer, const Self input, Index num_coeffs, half* output) {
237 const Index thread_id = blockIdx.x * blockDim.x + threadIdx.x;
238 const Index num_threads = blockDim.x * gridDim.x;
239 typedef typename packet_traits<Eigen::half>::type PacketType;
240 EIGEN_UNUSED_VARIABLE(input);
241
242 const Index num_packets = num_coeffs / Index(unpacket_traits<PacketType>::size);
243 PacketType* p_output = reinterpret_cast<PacketType*>(output);
244 for (Index i = thread_id; i < num_packets; i += num_threads) {
245 p_output[i] = reducer.template initializePacket<PacketType>();
246 }
247 Index packet_remainder = num_coeffs % Index(unpacket_traits<PacketType>::size);
248 if (thread_id < packet_remainder) {
249 output[num_coeffs - packet_remainder + thread_id] = reducer.initialize();
250 }
251}
252
253template <int BlockSize, int NumPerThread, typename Self,
254 typename Reducer, typename Index>
255__global__ void FullReductionKernelHalfFloat(Reducer reducer, const Self input, Index num_coeffs,
256 half* output, half2* scratch) {
257#if defined(__CUDA_ARCH__) && __CUDA_ARCH__ >= 300
258 EIGEN_UNUSED_VARIABLE(num_coeffs);
259 typedef typename packet_traits<Eigen::half>::type PacketType;
260 const int packet_width = unpacket_traits<PacketType>::size;
261 eigen_assert(NumPerThread % packet_width == 0);
262 const Index first_index = blockIdx.x * BlockSize * NumPerThread + packet_width * threadIdx.x;
263
264 // Initialize the output value if it wasn't initialized by the ReductionInitKernel
265 if (gridDim.x == 1) {
266 if (first_index == 0) {
267 int rem = num_coeffs % packet_width;
268 if (rem != 0) {
269 half2* p_scratch = reinterpret_cast<half2*>(scratch);
270 pstoreu(scratch, reducer.template initializePacket<PacketType>());
271 for (int i = 0; i < rem / 2; i++) {
272 *p_scratch = __halves2half2(input.coeff(num_coeffs - packet_width + 2 * i),
273 input.coeff(num_coeffs - packet_width + 2 * i + 1));
274 p_scratch++;
275 }
276 if ((num_coeffs & 1) != 0) {
277 half last = input.coeff(num_coeffs - 1);
278 *p_scratch = __halves2half2(last, reducer.initialize());
279 }
280 } else {
281 PacketType reduce = reducer.template initializePacket<PacketType>();
282 pstoreu(scratch, reduce);
283 }
284 }
285 __syncthreads();
286 }
287
288 PacketType accum = reducer.template initializePacket<PacketType>();
289 const Index max_iter =
290 numext::mini<Index>((num_coeffs - first_index) / packet_width, NumPerThread * BlockSize / packet_width);
291 for (Index i = 0; i < max_iter; i += BlockSize) {
292 const Index index = first_index + packet_width * i;
293 eigen_assert(index + packet_width < num_coeffs);
294 PacketType val = input.template packet<Unaligned>(index);
295 reducer.reducePacket(val, &accum);
296 }
297
298 #pragma unroll
299 for (int offset = warpSize / 2; offset > 0; offset /= 2) {
300#if defined(EIGEN_CUDA_SDK_VER) && EIGEN_CUDA_SDK_VER < 90000
301 PacketType r1;
302 half2* hr = reinterpret_cast<half2*>(&r1);
303 half2* hacc = reinterpret_cast<half2*>(&accum);
304 for (int i = 0; i < packet_width / 2; i++) {
305 hr[i] = __shfl_down(hacc[i], offset, warpSize);
306 }
307 reducer.reducePacket(r1, &accum);
308#else
309 PacketType r1;
310 half2* hr = reinterpret_cast<half2*>(&r1);
311 half2* hacc = reinterpret_cast<half2*>(&accum);
312 for (int i = 0; i < packet_width / 2; i++) {
313 hr[i] = __shfl_down_sync(0xFFFFFFFF, hacc[i], (unsigned)offset, warpSize);
314 }
315 reducer.reducePacket(r1, &accum);
316
317#endif
318 }
319
320 if ((threadIdx.x & (warpSize - 1)) == 0) {
321 atomicReduce(reinterpret_cast<PacketType*>(scratch), accum, reducer);
322 }
323
324 __syncthreads();
325 half2* rv1 = reinterpret_cast<half2*>(scratch);
326 if (packet_width > 2) {
327 reducer.reducePacket(rv1[2], rv1);
328 reducer.reducePacket(rv1[3], rv1 + 1);
329 reducer.reducePacket(rv1[1], rv1);
330 }
331 if (gridDim.x == 1) {
332 if (first_index == 0) {
333 half tmp = __low2half(*rv1);
334 reducer.reduce(__high2half(*rv1), &tmp);
335 *output = tmp;
336 }
337 }
338#else
339 EIGEN_UNUSED_VARIABLE(reducer);
340 EIGEN_UNUSED_VARIABLE(input);
341 EIGEN_UNUSED_VARIABLE(num_coeffs);
342 EIGEN_UNUSED_VARIABLE(output);
343 EIGEN_UNUSED_VARIABLE(scratch);
344#endif
345}
346
347template <typename Op>
348__global__ void ReductionCleanupKernelHalfFloat(Op& reducer, half* output, half2* scratch) {
349 eigen_assert(threadIdx.x == 1);
350 half tmp = __low2half(*scratch);
351 reducer.reduce(__high2half(*scratch), &tmp);
352 *output = tmp;
353}
354
355#endif
356
357template <typename Self, typename Op, typename OutputType, bool PacketAccess, typename Enabled = void>
358struct FullReductionLauncher {
359 static void run(const Self&, Op&, const GpuDevice&, OutputType*, typename Self::Index) {
360 assert(false && "Should only be called on doubles, floats and half floats");
361 }
362};
363
364// Specialization for float and double
365template <typename Self, typename Op, typename OutputType, bool PacketAccess>
366struct FullReductionLauncher<
367 Self, Op, OutputType, PacketAccess,
368 typename internal::enable_if<
369 internal::is_same<float, OutputType>::value ||
370 internal::is_same<double, OutputType>::value,
371 void>::type> {
372 static void run(const Self& self, Op& reducer, const GpuDevice& device, OutputType* output, typename Self::Index num_coeffs) {
373 typedef typename Self::Index Index;
374 const int block_size = 256;
375 const int num_per_thread = 128;
376 const int num_blocks = divup<int>(num_coeffs, block_size * num_per_thread);
377
378 unsigned int* semaphore = NULL;
379 if (num_blocks > 1) {
380 semaphore = device.semaphore();
381 }
382
383 LAUNCH_CUDA_KERNEL((FullReductionKernel<block_size, num_per_thread, Self, Op, Index>),
384 num_blocks, block_size, 0, device, reducer, self, num_coeffs, output, semaphore);
385 }
386};
387
388#ifdef EIGEN_HAS_CUDA_FP16
389template <typename Self, typename Op>
390struct FullReductionLauncher<Self, Op, Eigen::half, false> {
391 static void run(const Self&, Op&, const GpuDevice&, half*, typename Self::Index) {
392 assert(false && "Should not be called since there is no packet accessor");
393 }
394};
395
396template <typename Self, typename Op>
397struct FullReductionLauncher<Self, Op, Eigen::half, true> {
398 static void run(const Self& self, Op& reducer, const GpuDevice& device, half* output, typename Self::Index num_coeffs) {
399 typedef typename Self::Index Index;
400
401 const int block_size = 256;
402 const int num_per_thread = 128;
403 const int num_blocks = divup<int>(num_coeffs, block_size * num_per_thread);
404 half2* scratch = static_cast<half2*>(device.scratchpad());
405
406 if (num_blocks > 1) {
407 // We initialize the output and the scrathpad outside the reduction kernel when we can't be sure that there
408 // won't be a race conditions between multiple thread blocks.
409 LAUNCH_CUDA_KERNEL((ReductionInitFullReduxKernelHalfFloat<Self, Op, Index>),
410 1, 1, 0, device, reducer, self, num_coeffs, scratch);
411 }
412
413 LAUNCH_CUDA_KERNEL((FullReductionKernelHalfFloat<block_size, num_per_thread, Self, Op, Index>),
414 num_blocks, block_size, 0, device, reducer, self, num_coeffs, output, scratch);
415
416 if (num_blocks > 1) {
417 LAUNCH_CUDA_KERNEL((ReductionCleanupKernelHalfFloat<Op>),
418 1, 1, 0, device, reducer, output, scratch);
419 }
420 }
421};
422#endif
423
424
425template <typename Self, typename Op, bool Vectorizable>
426struct FullReducer<Self, Op, GpuDevice, Vectorizable> {
427 // Unfortunately nvidia doesn't support well exotic types such as complex,
428 // so reduce the scope of the optimized version of the code to the simple cases
429 // of doubles, floats and half floats
430#ifdef EIGEN_HAS_CUDA_FP16
431 static const bool HasOptimizedImplementation = !Op::IsStateful &&
432 (internal::is_same<typename Self::CoeffReturnType, float>::value ||
433 internal::is_same<typename Self::CoeffReturnType, double>::value ||
434 (internal::is_same<typename Self::CoeffReturnType, Eigen::half>::value && reducer_traits<Op, GpuDevice>::PacketAccess));
435#else
436 static const bool HasOptimizedImplementation = !Op::IsStateful &&
437 (internal::is_same<typename Self::CoeffReturnType, float>::value ||
438 internal::is_same<typename Self::CoeffReturnType, double>::value);
439#endif
440
441 template <typename OutputType>
442 static void run(const Self& self, Op& reducer, const GpuDevice& device, OutputType* output) {
443 assert(HasOptimizedImplementation && "Should only be called on doubles, floats or half floats");
444 const Index num_coeffs = array_prod(self.m_impl.dimensions());
445 // Don't crash when we're called with an input tensor of size 0.
446 if (num_coeffs == 0) {
447 return;
448 }
449
450 FullReductionLauncher<Self, Op, OutputType, reducer_traits<Op, GpuDevice>::PacketAccess>::run(self, reducer, device, output, num_coeffs);
451 }
452};
453
454
455template <int NumPerThread, typename Self,
456 typename Reducer, typename Index>
457__global__ void InnerReductionKernel(Reducer reducer, const Self input, Index num_coeffs_to_reduce, Index num_preserved_coeffs,
458 typename Self::CoeffReturnType* output) {
459#if __CUDA_ARCH__ >= 300
460 typedef typename Self::CoeffReturnType Type;
461 eigen_assert(blockDim.y == 1);
462 eigen_assert(blockDim.z == 1);
463 eigen_assert(gridDim.y == 1);
464 eigen_assert(gridDim.z == 1);
465
466 const int unroll_times = 16;
467 eigen_assert(NumPerThread % unroll_times == 0);
468
469 const Index input_col_blocks = divup<Index>(num_coeffs_to_reduce, blockDim.x * NumPerThread);
470 const Index num_input_blocks = input_col_blocks * num_preserved_coeffs;
471
472 const Index num_threads = blockDim.x * gridDim.x;
473 const Index thread_id = blockIdx.x * blockDim.x + threadIdx.x;
474
475 // Initialize the output values if they weren't initialized by the ReductionInitKernel
476 if (gridDim.x == 1) {
477 for (Index i = thread_id; i < num_preserved_coeffs; i += num_threads) {
478 output[i] = reducer.initialize();
479 }
480 __syncthreads();
481 }
482
483 for (Index i = blockIdx.x; i < num_input_blocks; i += gridDim.x) {
484 const Index row = i / input_col_blocks;
485
486 if (row < num_preserved_coeffs) {
487 const Index col_block = i % input_col_blocks;
488 const Index col_begin = col_block * blockDim.x * NumPerThread + threadIdx.x;
489
490 Type reduced_val = reducer.initialize();
491
492 for (Index j = 0; j < NumPerThread; j += unroll_times) {
493 const Index last_col = col_begin + blockDim.x * (j + unroll_times - 1);
494 if (last_col >= num_coeffs_to_reduce) {
495 for (Index col = col_begin + blockDim.x * j; col < num_coeffs_to_reduce; col += blockDim.x) {
496 const Type val = input.m_impl.coeff(row * num_coeffs_to_reduce + col);
497 reducer.reduce(val, &reduced_val);
498 }
499 break;
500 } else {
501 // Faster version of the loop with no branches after unrolling.
502#pragma unroll
503 for (int k = 0; k < unroll_times; ++k) {
504 const Index col = col_begin + blockDim.x * (j + k);
505 reducer.reduce(input.m_impl.coeff(row * num_coeffs_to_reduce + col), &reduced_val);
506 }
507 }
508 }
509
510#pragma unroll
511 for (int offset = warpSize/2; offset > 0; offset /= 2) {
512#if defined(EIGEN_CUDA_SDK_VER) && EIGEN_CUDA_SDK_VER < 90000
513 reducer.reduce(__shfl_down(reduced_val, offset), &reduced_val);
514#else
515 reducer.reduce(__shfl_down_sync(0xFFFFFFFF, reduced_val, offset), &reduced_val);
516#endif
517 }
518
519 if ((threadIdx.x & (warpSize - 1)) == 0) {
520 atomicReduce(&(output[row]), reduced_val, reducer);
521 }
522 }
523 }
524#else
525 EIGEN_UNUSED_VARIABLE(reducer);
526 EIGEN_UNUSED_VARIABLE(input);
527 EIGEN_UNUSED_VARIABLE(num_coeffs_to_reduce);
528 EIGEN_UNUSED_VARIABLE(num_preserved_coeffs);
529 EIGEN_UNUSED_VARIABLE(output);
530 assert(0 && "Shouldn't be called on unsupported device");
531#endif
532}
533
534#ifdef EIGEN_HAS_CUDA_FP16
535
536template <int NumPerThread, typename Self,
537 typename Reducer, typename Index>
538__global__ void InnerReductionKernelHalfFloat(Reducer reducer, const Self input, Index num_coeffs_to_reduce, Index num_preserved_coeffs,
539 half* output) {
540#if defined(__CUDA_ARCH__) && __CUDA_ARCH__ >= 300
541 eigen_assert(blockDim.y == 1);
542 eigen_assert(blockDim.z == 1);
543 eigen_assert(gridDim.y == 1);
544 eigen_assert(gridDim.z == 1);
545
546 typedef typename packet_traits<Eigen::half>::type PacketType;
547 const int packet_width = unpacket_traits<PacketType>::size;
548 const int unroll_times = 16 / packet_width;
549 eigen_assert(NumPerThread % unroll_times == 0);
550 eigen_assert(unroll_times % 2 == 0);
551
552 const Index input_col_blocks = numext::div_ceil<Index>(num_coeffs_to_reduce, blockDim.x * NumPerThread * 2);
553 const Index num_input_blocks = numext::div_ceil<Index>(input_col_blocks * num_preserved_coeffs, 2);
554
555 const Index num_threads = blockDim.x * gridDim.x;
556 const Index thread_id = blockIdx.x * blockDim.x + threadIdx.x;
557
558 // Initialize the output values if they weren't initialized by the ReductionInitKernel
559 if (gridDim.x == 1) {
560 Index i = packet_width * thread_id;
561 for (; i + packet_width <= num_preserved_coeffs; i += packet_width * num_threads) {
562 PacketType* poutput = reinterpret_cast<PacketType*>(output + i);
563 *poutput = reducer.template initializePacket<PacketType>();
564 }
565 if (i < num_preserved_coeffs) {
566 output[i] = reducer.initialize();
567 }
568 __syncthreads();
569 }
570
571 for (Index i = blockIdx.x; i < num_input_blocks; i += gridDim.x) {
572 const Index row = 2 * (i / input_col_blocks); // everybody takes 2 rows
573
574 if (row + 1 < num_preserved_coeffs) {
575 const Index col_block = i % input_col_blocks;
576 const Index col_begin = packet_width * (col_block * blockDim.x * NumPerThread + threadIdx.x);
577
578 PacketType reduced_val1 = reducer.template initializePacket<PacketType>();
579 PacketType reduced_val2 = reducer.template initializePacket<PacketType>();
580
581 for (Index j = 0; j < NumPerThread; j += unroll_times) {
582 const Index last_col = col_begin + blockDim.x * (j + unroll_times - 1) * packet_width;
583 if (last_col >= num_coeffs_to_reduce) {
584 Index col = col_begin + blockDim.x * j;
585 for (; col + packet_width <= num_coeffs_to_reduce; col += blockDim.x) {
586 const PacketType val1 = input.m_impl.template packet<Unaligned>(row * num_coeffs_to_reduce + col);
587 reducer.reducePacket(val1, &reduced_val1);
588 const PacketType val2 = input.m_impl.template packet<Unaligned>((row + 1) * num_coeffs_to_reduce + col);
589 reducer.reducePacket(val2, &reduced_val2);
590 }
591 if (col < num_coeffs_to_reduce) {
592 PacketType r1 = reducer.template initializePacket<PacketType>();
593 PacketType r2 = reducer.template initializePacket<PacketType>();
594 half2* hr1 = reinterpret_cast<half2*>(&r1);
595 half2* hr2 = reinterpret_cast<half2*>(&r2);
596 while (col + 1 < num_coeffs_to_reduce) {
597 *hr1 = __halves2half2(input.m_impl.coeff(row * num_coeffs_to_reduce + col),
598 input.m_impl.coeff(row * num_coeffs_to_reduce + col + 1));
599 *hr2 = __halves2half2(input.m_impl.coeff((row + 1) * num_coeffs_to_reduce + col),
600 input.m_impl.coeff((row + 1) * num_coeffs_to_reduce + col + 1));
601 hr1++;
602 hr2++;
603 col += 2;
604 }
605 if (col < num_coeffs_to_reduce) {
606 // Peel;
607 const half last1 = input.m_impl.coeff(row * num_coeffs_to_reduce + col);
608 *hr1 = __halves2half2(last1, reducer.initialize());
609 const half last2 = input.m_impl.coeff((row + 1) * num_coeffs_to_reduce + col);
610 *hr2 = __halves2half2(last2, reducer.initialize());
611 }
612 reducer.reducePacket(r1, &reduced_val1);
613 reducer.reducePacket(r2, &reduced_val2);
614 }
615 break;
616 } else {
617 // Faster version of the loop with no branches after unrolling.
618 #pragma unroll
619 for (int k = 0; k < unroll_times; ++k) {
620 const Index col = col_begin + blockDim.x * (j + k) * packet_width;
621 reducer.reducePacket(input.m_impl.template packet<Unaligned>(row * num_coeffs_to_reduce + col),
622 &reduced_val1);
623 reducer.reducePacket(input.m_impl.template packet<Unaligned>((row + 1) * num_coeffs_to_reduce + col),
624 &reduced_val2);
625 }
626 }
627 }
628
629 #pragma unroll
630 for (int offset = warpSize / 2; offset > 0; offset /= 2) {
631#if defined(EIGEN_CUDA_SDK_VER) && EIGEN_CUDA_SDK_VER < 90000
632 PacketType r1;
633 PacketType r2;
634 half2* hr1 = reinterpret_cast<half2*>(&r1);
635 half2* hr2 = reinterpret_cast<half2*>(&r2);
636 half2* rv1 = reinterpret_cast<half2*>(&reduced_val1);
637 half2* rv2 = reinterpret_cast<half2*>(&reduced_val2);
638 for (int i = 0; i < packet_width / 2; i++) {
639 hr1[i] = __shfl_down(rv1[i], offset, warpSize);
640 hr2[i] = __shfl_down(rv2[i], offset, warpSize);
641 }
642 reducer.reducePacket(r1, &reduced_val1);
643 reducer.reducePacket(r2, &reduced_val2);
644#else
645 PacketType r1;
646 PacketType r2;
647 half2* hr1 = reinterpret_cast<half2*>(&r1);
648 half2* hr2 = reinterpret_cast<half2*>(&r2);
649 half2* rr1 = reinterpret_cast<half2*>(&reduced_val1);
650 half2* rr2 = reinterpret_cast<half2*>(&reduced_val2);
651 for (int j = 0; j < packet_width / 2; j++) {
652 hr1[j] = __shfl_down_sync(0xFFFFFFFF, rr1[j], (unsigned)offset, warpSize);
653 hr2[j] = __shfl_down_sync(0xFFFFFFFF, rr2[j], (unsigned)offset, warpSize);
654 }
655 reducer.reducePacket(r1, &reduced_val1);
656 reducer.reducePacket(r2, &reduced_val2);
657
658#endif
659 }
660 half2* rv1 = reinterpret_cast<half2*>(&reduced_val1);
661 half2* rv2 = reinterpret_cast<half2*>(&reduced_val2);
662 half2 val;
663 if (packet_width > 2) {
664 reducer.reducePacket(rv1[2], rv1);
665 reducer.reducePacket(rv1[3], rv1 + 1);
666 reducer.reducePacket(rv1[1], rv1);
667 reducer.reducePacket(rv2[2], rv2);
668 reducer.reducePacket(rv2[3], rv2 + 1);
669 reducer.reducePacket(rv2[1], rv2);
670 }
671 half val1 = __low2half(*rv1);
672 reducer.reduce(__high2half(*rv1), &val1);
673 half val2 = __low2half(*rv2);
674 reducer.reduce(__high2half(*rv2), &val2);
675 val = __halves2half2(val1, val2);
676 if ((threadIdx.x & (warpSize - 1)) == 0) {
677 half* loc = output + row;
678 atomicReduce(reinterpret_cast<half2*>(loc), val, reducer);
679 }
680 }
681 }
682#else
683 EIGEN_UNUSED_VARIABLE(reducer);
684 EIGEN_UNUSED_VARIABLE(input);
685 EIGEN_UNUSED_VARIABLE(num_coeffs_to_reduce);
686 EIGEN_UNUSED_VARIABLE(num_preserved_coeffs);
687 EIGEN_UNUSED_VARIABLE(output);
688#endif
689}
690
691#endif
692
693template <typename Self, typename Op, typename OutputType, bool PacketAccess, typename Enabled = void>
694struct InnerReductionLauncher {
695 static EIGEN_DEVICE_FUNC bool run(const Self&, Op&, const GpuDevice&, OutputType*, typename Self::Index, typename Self::Index) {
696 assert(false && "Should only be called to reduce doubles, floats and half floats on a gpu device");
697 return true;
698 }
699};
700
701// Specialization for float and double
702template <typename Self, typename Op, typename OutputType, bool PacketAccess>
703struct InnerReductionLauncher<
704 Self, Op, OutputType, PacketAccess,
705 typename internal::enable_if<
706 internal::is_same<float, OutputType>::value ||
707 internal::is_same<double, OutputType>::value,
708 void>::type> {
709 static bool run(const Self& self, Op& reducer, const GpuDevice& device, OutputType* output, typename Self::Index num_coeffs_to_reduce, typename Self::Index num_preserved_vals) {
710 typedef typename Self::Index Index;
711
712 const Index num_coeffs = num_coeffs_to_reduce * num_preserved_vals;
713 const int block_size = 256;
714 const int num_per_thread = 128;
715 const int dyn_blocks = divup<int>(num_coeffs, block_size * num_per_thread);
716 const int max_blocks = device.getNumCudaMultiProcessors() *
717 device.maxCudaThreadsPerMultiProcessor() / block_size;
718 const int num_blocks = numext::mini<int>(max_blocks, dyn_blocks);
719
720 if (num_blocks > 1) {
721 // We initialize the outputs outside the reduction kernel when we can't be sure that there
722 // won't be a race conditions between multiple thread blocks.
723 const int dyn_blocks2 = divup<int>(num_preserved_vals, 1024);
724 const int max_blocks2 = device.getNumCudaMultiProcessors() *
725 device.maxCudaThreadsPerMultiProcessor() / 1024;
726 const int num_blocks2 = numext::mini<int>(max_blocks2, dyn_blocks2);
727 LAUNCH_CUDA_KERNEL((ReductionInitKernel<OutputType, Index>),
728 num_blocks2, 1024, 0, device, reducer.initialize(),
729 num_preserved_vals, output);
730 }
731
732 LAUNCH_CUDA_KERNEL((InnerReductionKernel<num_per_thread, Self, Op, Index>),
733 num_blocks, block_size, 0, device, reducer, self, num_coeffs_to_reduce, num_preserved_vals, output);
734
735 return false;
736 }
737};
738
739#ifdef EIGEN_HAS_CUDA_FP16
740template <typename Self, typename Op>
741struct InnerReductionLauncher<Self, Op, Eigen::half, false> {
742 static bool run(const Self&, Op&, const GpuDevice&, half*, typename Self::Index, typename Self::Index) {
743 assert(false && "Should not be called since there is no packet accessor");
744 return true;
745 }
746};
747
748template <typename Self, typename Op>
749struct InnerReductionLauncher<Self, Op, Eigen::half, true> {
750 static bool run(const Self& self, Op& reducer, const GpuDevice& device, half* output, typename Self::Index num_coeffs_to_reduce, typename Self::Index num_preserved_vals) {
751 typedef typename Self::Index Index;
752
753 if (num_preserved_vals % 2 != 0) {
754 // Not supported yet, revert to the slower code path
755 return true;
756 }
757
758 const Index num_coeffs = num_coeffs_to_reduce * num_preserved_vals;
759 const int block_size = /*256*/128;
760 const int num_per_thread = /*128*/64;
761 const int dyn_blocks = divup<int>(num_coeffs, block_size * num_per_thread);
762 const int max_blocks = device.getNumCudaMultiProcessors() *
763 device.maxCudaThreadsPerMultiProcessor() / block_size;
764 const int num_blocks = numext::mini<int>(max_blocks, dyn_blocks);
765
766 if (num_blocks > 1) {
767 // We initialize the outputs outside the reduction kernel when we can't be sure that there
768 // won't be a race conditions between multiple thread blocks.
769 LAUNCH_CUDA_KERNEL((ReductionInitKernelHalfFloat<Self, Op, Index>),
770 1, 1, 0, device, reducer, self, num_preserved_vals, output);
771 }
772
773 LAUNCH_CUDA_KERNEL((InnerReductionKernelHalfFloat<num_per_thread, Self, Op, Index>),
774 num_blocks, block_size, 0, device, reducer, self, num_coeffs_to_reduce, num_preserved_vals, output);
775
776 return false;
777 }
778};
779#endif
780
781
782template <typename Self, typename Op>
783struct InnerReducer<Self, Op, GpuDevice> {
784 // Unfortunately nvidia doesn't support well exotic types such as complex,
785 // so reduce the scope of the optimized version of the code to the simple case
786 // of floats and half floats.
787#ifdef EIGEN_HAS_CUDA_FP16
788 static const bool HasOptimizedImplementation = !Op::IsStateful &&
789 (internal::is_same<typename Self::CoeffReturnType, float>::value ||
790 internal::is_same<typename Self::CoeffReturnType, double>::value ||
791 (internal::is_same<typename Self::CoeffReturnType, Eigen::half>::value && reducer_traits<Op, GpuDevice>::PacketAccess));
792#else
793 static const bool HasOptimizedImplementation = !Op::IsStateful &&
794 (internal::is_same<typename Self::CoeffReturnType, float>::value ||
795 internal::is_same<typename Self::CoeffReturnType, double>::value);
796#endif
797
798 template <typename OutputType>
799 static bool run(const Self& self, Op& reducer, const GpuDevice& device, OutputType* output, typename Self::Index num_coeffs_to_reduce, typename Self::Index num_preserved_vals) {
800 assert(HasOptimizedImplementation && "Should only be called on doubles, floats or half floats");
801 const Index num_coeffs = array_prod(self.m_impl.dimensions());
802 // Don't crash when we're called with an input tensor of size 0.
803 if (num_coeffs == 0) {
804 return true;
805 }
806 // It's faster to use the usual code.
807 if (num_coeffs_to_reduce <= 128) {
808 return true;
809 }
810
811 return InnerReductionLauncher<Self, Op, OutputType, reducer_traits<Op, GpuDevice>::PacketAccess>::run(self, reducer, device, output, num_coeffs_to_reduce, num_preserved_vals);
812 }
813};
814
815template <int NumPerThread, typename Self,
816 typename Reducer, typename Index>
817__global__ void OuterReductionKernel(Reducer reducer, const Self input, Index num_coeffs_to_reduce, Index num_preserved_coeffs,
818 typename Self::CoeffReturnType* output) {
819 const Index num_threads = blockDim.x * gridDim.x;
820 const Index thread_id = blockIdx.x * blockDim.x + threadIdx.x;
821 // Initialize the output values if they weren't initialized by the ReductionInitKernel
822 if (gridDim.x == 1) {
823 for (Index i = thread_id; i < num_preserved_coeffs; i += num_threads) {
824 output[i] = reducer.initialize();
825 }
826 __syncthreads();
827 }
828
829 // Do the reduction.
830 const Index max_iter = num_preserved_coeffs * divup<Index>(num_coeffs_to_reduce, NumPerThread);
831 for (Index i = thread_id; i < max_iter; i += num_threads) {
832 const Index input_col = i % num_preserved_coeffs;
833 const Index input_row = (i / num_preserved_coeffs) * NumPerThread;
834 typename Self::CoeffReturnType reduced_val = reducer.initialize();
835 const Index max_row = numext::mini(input_row + NumPerThread, num_coeffs_to_reduce);
836 for (Index j = input_row; j < max_row; j++) {
837 typename Self::CoeffReturnType val = input.m_impl.coeff(j * num_preserved_coeffs + input_col);
838 reducer.reduce(val, &reduced_val);
839 }
840 atomicReduce(&(output[input_col]), reduced_val, reducer);
841 }
842}
843
844
845template <typename Self, typename Op>
846struct OuterReducer<Self, Op, GpuDevice> {
847 // Unfortunately nvidia doesn't support well exotic types such as complex,
848 // so reduce the scope of the optimized version of the code to the simple case
849 // of floats.
850 static const bool HasOptimizedImplementation = !Op::IsStateful &&
851 (internal::is_same<typename Self::CoeffReturnType, float>::value ||
852 internal::is_same<typename Self::CoeffReturnType, double>::value);
853 template <typename Device, typename OutputType>
854 static EIGEN_DEVICE_FUNC bool run(const Self&, Op&, const Device&, OutputType*, typename Self::Index, typename Self::Index) {
855 assert(false && "Should only be called to reduce doubles or floats on a gpu device");
856 return true;
857 }
858
859 static bool run(const Self& self, Op& reducer, const GpuDevice& device, float* output, typename Self::Index num_coeffs_to_reduce, typename Self::Index num_preserved_vals) {
860 typedef typename Self::Index Index;
861
862 // It's faster to use the usual code.
863 if (num_coeffs_to_reduce <= 32) {
864 return true;
865 }
866
867 const Index num_coeffs = num_coeffs_to_reduce * num_preserved_vals;
868 const int block_size = 256;
869 const int num_per_thread = 16;
870 const int dyn_blocks = divup<int>(num_coeffs, block_size * num_per_thread);
871 const int max_blocks = device.getNumCudaMultiProcessors() *
872 device.maxCudaThreadsPerMultiProcessor() / block_size;
873 const int num_blocks = numext::mini<int>(max_blocks, dyn_blocks);
874
875 if (num_blocks > 1) {
876 // We initialize the outputs in the reduction kernel itself when we don't have to worry
877 // about race conditions between multiple thread blocks.
878 const int dyn_blocks2 = divup<int>(num_preserved_vals, 1024);
879 const int max_blocks2 = device.getNumCudaMultiProcessors() *
880 device.maxCudaThreadsPerMultiProcessor() / 1024;
881 const int num_blocks2 = numext::mini<int>(max_blocks2, dyn_blocks2);
882 LAUNCH_CUDA_KERNEL((ReductionInitKernel<float, Index>),
883 num_blocks2, 1024, 0, device, reducer.initialize(),
884 num_preserved_vals, output);
885 }
886
887 LAUNCH_CUDA_KERNEL((OuterReductionKernel<num_per_thread, Self, Op, Index>),
888 num_blocks, block_size, 0, device, reducer, self, num_coeffs_to_reduce, num_preserved_vals, output);
889
890 return false;
891 }
892};
893
894#endif
895
896
897} // end namespace internal
898} // end namespace Eigen
899
900#endif // EIGEN_CXX11_TENSOR_TENSOR_REDUCTION_CUDA_H
Namespace containing all symbols from the Eigen library.
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index