SkePU(integratedwithStarPU)  0.8.1
 All Classes Namespaces Files Functions Enumerations Friends Macros Groups Pages
map_kernels.h
Go to the documentation of this file.
1 
5 #ifndef MAP_KERNELS_H
6 #define MAP_KERNELS_H
7 
8 #ifdef SKEPU_OPENCL
9 
10 #include <string>
11 
12 namespace skepu
13 {
14 
30 static std::string UnaryMapKernel_CL(
31 "__kernel void UnaryMapKernel_KERNELNAME(__global TYPE* input, __global TYPE* output, unsigned int numElements, unsigned int offset1, unsigned int offset2, CONST_TYPE const1)\n"
32 "{\n"
33 " input = (__global void *)input + offset1; /* partitioning is special with opencl */ \n"
34 " output = (__global void *)output + offset2; /* partitioning is special with opencl */ \n"
35 " int i = get_global_id(0);\n"
36 " unsigned int gridSize = get_local_size(0)*get_num_groups(0);\n"
37 " while(i < numElements)\n"
38 " {\n"
39 " output[i] = FUNCTIONNAME(input[i], const1);\n"
40 " i += gridSize;\n"
41 " }\n"
42 "}\n"
43 );
44 
49 static std::string BinaryMapKernel_CL(
50 "__kernel void BinaryMapKernel_KERNELNAME(__global TYPE* input1, __global TYPE* input2, __global TYPE* output, unsigned int n, unsigned int offset1, unsigned int offset2, unsigned int offset3, CONST_TYPE const1)\n"
51 "{\n"
52 " input1 = (__global void *)input1 + offset1; /* partitioning is special with opencl */ \n"
53 " input2 = (__global void *)input2 + offset2; /* partitioning is special with opencl */ \n"
54 " output = (__global void *)output + offset3; /* partitioning is special with opencl */ \n"
55 " int i = get_global_id(0);\n"
56 " unsigned int gridSize = get_local_size(0)*get_num_groups(0);\n"
57 " while(i < n)\n"
58 " {\n"
59 " output[i] = FUNCTIONNAME(input1[i], input2[i], const1);\n"
60 " i += gridSize;\n"
61 " }\n"
62 "}\n"
63 );
64 
69 static std::string TrinaryMapKernel_CL(
70 "__kernel void TrinaryMapKernel_KERNELNAME(__global TYPE* input1, __global TYPE* input2, __global TYPE* input3, __global TYPE* output, unsigned int n, unsigned int offset1, unsigned int offset2, unsigned int offset3, unsigned int offset4, CONST_TYPE const1)\n"
71 "{\n"
72 " input1 = (__global void *)input1 + offset1; /* partitioning is special with opencl */ \n"
73 " input2 = (__global void *)input2 + offset2; /* partitioning is special with opencl */ \n"
74 " input3 = (__global void *)input3 + offset3; /* partitioning is special with opencl */ \n"
75 " output = (__global void *)output + offset4; /* partitioning is special with opencl */ \n"
76 " int i = get_global_id(0);\n"
77 " unsigned int gridSize = get_local_size(0)*get_num_groups(0);\n"
78 " while(i < n)\n"
79 " {\n"
80 " output[i] = FUNCTIONNAME(input1[i], input2[i], input3[i], const1);\n"
81 " i += gridSize;\n"
82 " }\n"
83 "}\n"
84 );
85 
90 }
91 
92 #endif
93 
94 #ifdef SKEPU_CUDA
95 
96 namespace skepu
97 {
98 
108 template <typename T, typename U, typename UnaryFunc>
109 __global__ void MapKernelUnary_CU(UnaryFunc mapFunc, T* input, T* output, unsigned int n, U const1)
110 {
111  unsigned int i = blockIdx.x * blockDim.x + threadIdx.x;
112  unsigned int gridSize = blockDim.x*gridDim.x;
113 
114  while(i < n)
115  {
116  output[i] = mapFunc.CU(input[i], const1);
117  i += gridSize;
118  }
119 }
120 
125 template <typename T, typename U, typename BinaryFunc>
126 __global__ void MapKernelBinary_CU(BinaryFunc mapFunc, T* input1, T* input2, T* output, unsigned int n, U const1)
127 {
128  unsigned int i = blockIdx.x * blockDim.x + threadIdx.x;
129  unsigned int gridSize = blockDim.x*gridDim.x;
130 
131  while(i < n)
132  {
133  output[i] = mapFunc.CU(input1[i], input2[i], const1);
134  i += gridSize;
135  }
136 }
137 
142 template <typename T, typename U, typename TrinaryFunc>
143 __global__ void MapKernelTrinary_CU(TrinaryFunc mapFunc, T* input1, T* input2, T* input3, T* output, unsigned int n, U const1)
144 {
145  unsigned int i = blockIdx.x * blockDim.x + threadIdx.x;
146  unsigned int gridSize = blockDim.x*gridDim.x;
147 
148  while(i < n)
149  {
150  output[i] = mapFunc.CU(input1[i], input2[i], input3[i], const1);
151  i += gridSize;
152  }
153 }
154 
159 }
160 
161 #endif
162 
163 #endif
164 
__global__ void MapKernelUnary_CU(UnaryFunc mapFunc, T *input, T *output, unsigned int n, U const1)
Definition: map_kernels.h:109
static std::string UnaryMapKernel_CL("__kernel void UnaryMapKernel_KERNELNAME(__global TYPE* input, __global TYPE* output, unsigned int numElements, unsigned int offset1, unsigned int offset2, CONST_TYPE const1)\n""{\n"" input = (__global void *)input + offset1; /* partitioning is special with opencl */ \n"" output = (__global void *)output + offset2; /* partitioning is special with opencl */ \n"" int i = get_global_id(0);\n"" unsigned int gridSize = get_local_size(0)*get_num_groups(0);\n"" while(i < numElements)\n"" {\n"" output[i] = FUNCTIONNAME(input[i], const1);\n"" i += gridSize;\n"" }\n""}\n")
__global__ void MapKernelBinary_CU(BinaryFunc mapFunc, T *input1, T *input2, T *output, unsigned int n, U const1)
Definition: map_kernels.h:126
__global__ void MapKernelTrinary_CU(TrinaryFunc mapFunc, T *input1, T *input2, T *input3, T *output, unsigned int n, U const1)
Definition: map_kernels.h:143
static std::string BinaryMapKernel_CL("__kernel void BinaryMapKernel_KERNELNAME(__global TYPE* input1, __global TYPE* input2, __global TYPE* output, unsigned int n, unsigned int offset1, unsigned int offset2, unsigned int offset3, CONST_TYPE const1)\n""{\n"" input1 = (__global void *)input1 + offset1; /* partitioning is special with opencl */ \n"" input2 = (__global void *)input2 + offset2; /* partitioning is special with opencl */ \n"" output = (__global void *)output + offset3; /* partitioning is special with opencl */ \n"" int i = get_global_id(0);\n"" unsigned int gridSize = get_local_size(0)*get_num_groups(0);\n"" while(i < n)\n"" {\n"" output[i] = FUNCTIONNAME(input1[i], input2[i], const1);\n"" i += gridSize;\n"" }\n""}\n")
static std::string TrinaryMapKernel_CL("__kernel void TrinaryMapKernel_KERNELNAME(__global TYPE* input1, __global TYPE* input2, __global TYPE* input3, __global TYPE* output, unsigned int n, unsigned int offset1, unsigned int offset2, unsigned int offset3, unsigned int offset4, CONST_TYPE const1)\n""{\n"" input1 = (__global void *)input1 + offset1; /* partitioning is special with opencl */ \n"" input2 = (__global void *)input2 + offset2; /* partitioning is special with opencl */ \n"" input3 = (__global void *)input3 + offset3; /* partitioning is special with opencl */ \n"" output = (__global void *)output + offset4; /* partitioning is special with opencl */ \n"" int i = get_global_id(0);\n"" unsigned int gridSize = get_local_size(0)*get_num_groups(0);\n"" while(i < n)\n"" {\n"" output[i] = FUNCTIONNAME(input1[i], input2[i], input3[i], const1);\n"" i += gridSize;\n"" }\n""}\n")