|
SkePU 0.7
|
00001 00005 #ifndef MAP_KERNELS_H 00006 #define MAP_KERNELS_H 00007 00008 #ifdef SKEPU_OPENCL 00009 00010 #include <string> 00011 00012 namespace skepu 00013 { 00014 00031 static std::string UnaryMapKernel_CL( 00032 "__kernel void UnaryMapKernel_KERNELNAME(__global TYPE* input, __global TYPE* output, unsigned int numElements, TYPE const1)\n" 00033 "{\n" 00034 " int i = get_global_id(0);\n" 00035 " unsigned int gridSize = get_local_size(0)*get_num_groups(0);\n" 00036 " while(i < numElements)\n" 00037 " {\n" 00038 " output[i] = FUNCTIONNAME(input[i], const1);\n" 00039 " i += gridSize;\n" 00040 " }\n" 00041 "}\n" 00042 ); 00043 00049 static std::string BinaryMapKernel_CL( 00050 "__kernel void BinaryMapKernel_KERNELNAME(__global TYPE* input1, __global TYPE* input2, __global TYPE* output, unsigned int n, TYPE const1)\n" 00051 "{\n" 00052 " int i = get_global_id(0);\n" 00053 " unsigned int gridSize = get_local_size(0)*get_num_groups(0);\n" 00054 " while(i < n)\n" 00055 " {\n" 00056 " output[i] = FUNCTIONNAME(input1[i], input2[i], const1);\n" 00057 " i += gridSize;\n" 00058 " }\n" 00059 "}\n" 00060 ); 00061 00067 static std::string TrinaryMapKernel_CL( 00068 "__kernel void TrinaryMapKernel_KERNELNAME(__global TYPE* input1, __global TYPE* input2, __global TYPE* input3, __global TYPE* output, unsigned int n, TYPE const1)\n" 00069 "{\n" 00070 " int i = get_global_id(0);\n" 00071 " unsigned int gridSize = get_local_size(0)*get_num_groups(0);\n" 00072 " while(i < n)\n" 00073 " {\n" 00074 " output[i] = FUNCTIONNAME(input1[i], input2[i], input3[i], const1);\n" 00075 " i += gridSize;\n" 00076 " }\n" 00077 "}\n" 00078 ); 00079 00084 } 00085 00086 #endif 00087 00088 #ifdef SKEPU_CUDA 00089 00090 namespace skepu 00091 { 00092 00109 template <typename T, typename UnaryFunc> 00110 __global__ void MapKernelUnary_CU(UnaryFunc mapFunc, T* input, T* output, unsigned int n, T const1) 00111 { 00112 unsigned int i = blockIdx.x * blockDim.x + threadIdx.x; 00113 unsigned int gridSize = blockDim.x*gridDim.x; 00114 00115 while(i < n) 00116 { 00117 output[i] = mapFunc.CU(input[i], const1); 00118 i += gridSize; 00119 } 00120 } 00121 00127 template <typename T, typename BinaryFunc> 00128 __global__ void MapKernelBinary_CU(BinaryFunc mapFunc, T* input1, T* input2, T* output, unsigned int n, T const1) 00129 { 00130 unsigned int i = blockIdx.x * blockDim.x + threadIdx.x; 00131 unsigned int gridSize = blockDim.x*gridDim.x; 00132 00133 while(i < n) 00134 { 00135 output[i] = mapFunc.CU(input1[i], input2[i], const1); 00136 i += gridSize; 00137 } 00138 } 00139 00145 template <typename T, typename TrinaryFunc> 00146 __global__ void MapKernelTrinary_CU(TrinaryFunc mapFunc, T* input1, T* input2, T* input3, T* output, unsigned int n, T const1) 00147 { 00148 unsigned int i = blockIdx.x * blockDim.x + threadIdx.x; 00149 unsigned int gridSize = blockDim.x*gridDim.x; 00150 00151 while(i < n) 00152 { 00153 output[i] = mapFunc.CU(input1[i], input2[i], input3[i], const1); 00154 i += gridSize; 00155 } 00156 } 00157 00162 } 00163 00164 #endif 00165 00166 #endif 00167
1.7.4