ViennaCL - The Vienna Computing Library  1.7.1
Free open-source GPU-accelerated linear algebra and solver library.
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
custom-kernels.cpp
Go to the documentation of this file.
1 /* =========================================================================
2  Copyright (c) 2010-2016, Institute for Microelectronics,
3  Institute for Analysis and Scientific Computing,
4  TU Wien.
5  Portions of this software are copyright by UChicago Argonne, LLC.
6 
7  -----------------
8  ViennaCL - The Vienna Computing Library
9  -----------------
10 
11  Project Head: Karl Rupp rupp@iue.tuwien.ac.at
12 
13  (A list of authors and contributors can be found in the PDF manual)
14 
15  License: MIT (X11), see file LICENSE in the base directory
16 ============================================================================= */
17 
26 // System headers
27 #include <iostream>
28 #include <string>
29 
30 #ifndef VIENNACL_WITH_OPENCL
31  #define VIENNACL_WITH_OPENCL
32 #endif
33 
34 // ViennaCL headers
35 #include "viennacl/ocl/backend.hpp"
36 #include "viennacl/vector.hpp"
38 
39 
54 static const char * my_compute_program =
55 "__kernel void elementwise_prod(\n"
56 " __global const float * vec1,\n"
57 " __global const float * vec2, \n"
58 " __global float * result,\n"
59 " unsigned int size) \n"
60 "{ \n"
61 " for (unsigned int i = get_global_id(0); i < size; i += get_global_size(0))\n"
62 " result[i] = vec1[i] * vec2[i];\n"
63 "};\n\n"
64 "__kernel void elementwise_div(\n"
65 " __global const float * vec1,\n"
66 " __global const float * vec2, \n"
67 " __global float * result,\n"
68 " unsigned int size) \n"
69 "{ \n"
70 " for (unsigned int i = get_global_id(0); i < size; i += get_global_size(0))\n"
71 " result[i] = vec1[i] / vec2[i];\n"
72 "};\n";
73 
77 int main()
78 {
79  typedef float ScalarType;
80 
84  unsigned int vector_size = 10;
85  viennacl::vector<ScalarType> vec1(vector_size);
86  viennacl::vector<ScalarType> vec2(vector_size);
87  viennacl::vector<ScalarType> result_mul(vector_size);
88  viennacl::vector<ScalarType> result_div(vector_size);
89 
93  for (unsigned int i=0; i<vector_size; ++i)
94  {
95  vec1[i] = static_cast<ScalarType>(i);
96  vec2[i] = static_cast<ScalarType>(vector_size-i);
97  }
98 
103  viennacl::ocl::program & my_prog = viennacl::ocl::current_context().add_program(my_compute_program, "my_compute_program");
104  // Note: Releases older than ViennaCL 1.5.0 required calls to add_kernel(). This is no longer needed, the respective interface has been removed.
105 
111  viennacl::ocl::kernel & my_kernel_mul = my_prog.get_kernel("elementwise_prod");
112  viennacl::ocl::kernel & my_kernel_div = my_prog.get_kernel("elementwise_div");
113 
118  viennacl::ocl::enqueue(my_kernel_mul(vec1, vec2, result_mul, static_cast<cl_uint>(vec1.size())));
119  viennacl::ocl::enqueue(my_kernel_div(vec1, vec2, result_div, static_cast<cl_uint>(vec1.size())));
120 
124  std::cout << " vec1: " << vec1 << std::endl;
125  std::cout << " vec2: " << vec2 << std::endl;
126  std::cout << "vec1 .* vec2: " << result_mul << std::endl;
127  std::cout << "vec1 /* vec2: " << result_div << std::endl;
128  std::cout << "norm_2(vec1 .* vec2): " << viennacl::linalg::norm_2(result_mul) << std::endl;
129  std::cout << "norm_2(vec1 /* vec2): " << viennacl::linalg::norm_2(result_div) << std::endl;
130 
134  std::cout << "!!!! TUTORIAL COMPLETED SUCCESSFULLY !!!!" << std::endl;
135 
136  return EXIT_SUCCESS;
137 }
138 
T norm_2(std::vector< T, A > const &v1)
Definition: norm_2.hpp:96
Generic interface for the l^2-norm. See viennacl/linalg/vector_operations.hpp for implementations...
Represents an OpenCL kernel within ViennaCL.
Definition: kernel.hpp:58
int main()
Definition: bisect.cpp:91
viennacl::ocl::context & current_context()
Convenience function for returning the current context.
Definition: backend.hpp:213
viennacl::ocl::program & add_program(cl_program p, std::string const &prog_name)
Adds a program to the context.
Definition: context.hpp:368
Wrapper class for an OpenCL program.
Definition: program.hpp:42
Implementations of the OpenCL backend, where all contexts are stored in.
void enqueue(KernelType &k, viennacl::ocl::command_queue const &queue)
Enqueues a kernel in the provided queue.
Definition: enqueue.hpp:50
The vector type with operator-overloads and proxy classes is defined here. Linear algebra operations ...
float ScalarType
Definition: fft_1d.cpp:42
viennacl::ocl::kernel & get_kernel(std::string const &name)
Returns the kernel with the provided name.
Definition: context.hpp:773