File: cudax_kernels.cu

package info (click to toggle)
starpu-contrib 1.0.1%2Bdfsg-1
  • links: PTS, VCS
  • area: contrib
  • in suites: wheezy
  • size: 13,836 kB
  • sloc: ansic: 77,357; cpp: 23,334; sh: 12,088; makefile: 2,086; lisp: 758; yacc: 185; sed: 126; fortran: 13
file content (156 lines) | stat: -rw-r--r-- 4,929 bytes parent folder | download | duplicates (4)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
/* StarPU --- Runtime system for heterogeneous multicore architectures.
 *
 * Copyright (C) 2009, 2010  Université de Bordeaux 1
 * Copyright (C) 2010, 2011  Centre National de la Recherche Scientifique
 *
 * StarPU is free software; you can redistribute it and/or modify
 * it under the terms of the GNU Lesser General Public License as published by
 * the Free Software Foundation; either version 2.1 of the License, or (at
 * your option) any later version.
 *
 * StarPU is distributed in the hope that it will be useful, but
 * WITHOUT ANY WARRANTY; without even the implied warranty of
 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
 *
 * See the GNU Lesser General Public License in COPYING.LGPL for more details.
 */

#define _externC extern "C"
#include "cudax_kernels.h"

/* Note: these assume that the sizes are powers of two */

#define VARS_1d \
	unsigned start = threadIdx.x + blockIdx.x * blockDim.x; \
	unsigned numthreads = blockDim.x * gridDim.x;

#define DISTRIB_1d(n, func,args) \
	unsigned threads_per_block = 128; \
\
	if (n < threads_per_block) \
	{			   \
		dim3 dimGrid(n); \
		func <<<dimGrid, 1, 0, starpu_cuda_get_local_stream()>>> args; \
	} 					\
	else 					\
	{				     \
		dim3 dimGrid(n / threads_per_block); \
		dim3 dimBlock(threads_per_block); \
		func <<<dimGrid, dimBlock, 0, starpu_cuda_get_local_stream()>>> args; \
	} \
	cudaStreamSynchronize(starpu_cuda_get_local_stream()); \

extern "C" __global__ void
STARPUFFT(cuda_twist1_1d)(const _cuComplex *in, _cuComplex *twisted1, unsigned i, unsigned n1, unsigned n2)
{
	unsigned j;
	VARS_1d
	unsigned end = n2;

	for (j = start; j < end; j += numthreads)
		twisted1[j] = in[i+j*n1];
}

extern "C" void
STARPUFFT(cuda_twist1_1d_host)(const _cuComplex *in, _cuComplex *twisted1, unsigned i, unsigned n1, unsigned n2)
{
	DISTRIB_1d(n2, STARPUFFT(cuda_twist1_1d), (in, twisted1, i, n1, n2));
}

extern "C" __global__ void
STARPUFFT(cuda_twiddle_1d)(_cuComplex * out, const _cuComplex * roots, unsigned n, unsigned i)
{
	unsigned j;
	VARS_1d
	unsigned end = n;

	for (j = start; j < end; j += numthreads)
		out[j] = _cuCmul(out[j], roots[i*j]);
	return;
}

extern "C" void
STARPUFFT(cuda_twiddle_1d_host)(_cuComplex *out, const _cuComplex *roots, unsigned n, unsigned i)
{
	DISTRIB_1d(n, STARPUFFT(cuda_twiddle_1d), (out, roots, n, i));
}

#define VARS_2d \
	unsigned startx = threadIdx.x + blockIdx.x * blockDim.x; \
	unsigned starty = threadIdx.y + blockIdx.y * blockDim.y; \
	unsigned numthreadsx = blockDim.x * gridDim.x; \
	unsigned numthreadsy = blockDim.y * gridDim.y;

/* FIXME: introduce threads_per_dim_n / m instead */
#define DISTRIB_2d(n, m, func, args) \
	unsigned threads_per_dim = 16; \
	if (n < threads_per_dim) \
	{				   \
		if (m < threads_per_dim) \
		{			    \
			dim3 dimGrid(n, m); \
			func <<<dimGrid, 1, 0, starpu_cuda_get_local_stream()>>> args; \
		} \
		else \
		{					      \
			dim3 dimGrid(1, m / threads_per_dim); \
			dim3 dimBlock(n, threads_per_dim); \
			func <<<dimGrid, dimBlock, 0, starpu_cuda_get_local_stream()>>> args; \
		} \
	} \
	else \
	{				   \
		if (m < threads_per_dim) \
		{					      \
			dim3 dimGrid(n / threads_per_dim, 1); \
			dim3 dimBlock(threads_per_dim, m); \
			func <<<dimGrid, dimBlock, 0, starpu_cuda_get_local_stream()>>> args; \
		} \
		else \
		{							\
			dim3 dimGrid(n / threads_per_dim, m / threads_per_dim); \
			dim3 dimBlock(threads_per_dim, threads_per_dim); \
			func <<<dimGrid, dimBlock, 0, starpu_cuda_get_local_stream()>>> args; \
		} \
	} \
	cudaStreamSynchronize(starpu_cuda_get_local_stream()); \

extern "C" __global__ void
STARPUFFT(cuda_twist1_2d)(const _cuComplex *in, _cuComplex *twisted1, unsigned i, unsigned j, unsigned n1, unsigned n2, unsigned m1, unsigned m2)
{
	unsigned k, l;
	VARS_2d
	unsigned endx = n2;
	unsigned endy = m2;
	unsigned m = m1*m2;

	for (k = startx; k < endx; k += numthreadsx)
		for (l = starty; l < endy; l += numthreadsy)
			twisted1[k*m2+l] = in[i*m+j+k*m*n1+l*m1];
}

extern "C" void
STARPUFFT(cuda_twist1_2d_host)(const _cuComplex *in, _cuComplex *twisted1, unsigned i, unsigned j, unsigned n1, unsigned n2, unsigned m1, unsigned m2)
{
	DISTRIB_2d(n2, m2, STARPUFFT(cuda_twist1_2d), (in, twisted1, i, j, n1, n2, m1, m2));
}

extern "C" __global__ void
STARPUFFT(cuda_twiddle_2d)(_cuComplex * out, const _cuComplex * roots0, const _cuComplex * roots1, unsigned n2, unsigned m2, unsigned i, unsigned j)
{
	unsigned k, l;
	VARS_2d
	unsigned endx = n2;
	unsigned endy = m2;

	for (k = startx; k < endx ; k += numthreadsx)
		for (l = starty; l < endy ; l += numthreadsy)
			out[k*m2 + l] = _cuCmul(_cuCmul(out[k*m2 + l], roots0[i*k]), roots1[j*l]);
	return;
}

extern "C" void
STARPUFFT(cuda_twiddle_2d_host)(_cuComplex *out, const _cuComplex *roots0, const _cuComplex *roots1, unsigned n2, unsigned m2, unsigned i, unsigned j)
{
	DISTRIB_2d(n2, m2, STARPUFFT(cuda_twiddle_2d), (out, roots0, roots1, n2, m2, i, j));
}