ViennaCL - The Vienna Computing Library  1.7.0
Free open-source GPU-accelerated linear algebra and solver library.
libviennacl.cpp
Go to the documentation of this file.
1 /* =========================================================================
2  Copyright (c) 2010-2015, 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 
27 // include necessary system headers
28 #include <iostream>
29 #include <vector>
30 
31 // Some helper functions for this tutorial:
32 #include "viennacl.hpp"
33 
34 #include "viennacl/vector.hpp"
35 
36 
40 int main()
41 {
42  std::size_t size = 10;
43 
44  ViennaCLInt half_size = static_cast<ViennaCLInt>(size / 2);
45 
46 
51  ViennaCLBackend my_backend;
52  ViennaCLBackendCreate(&my_backend);
53 
54 
63 
64  ViennaCLHostSswap(my_backend, half_size,
65  viennacl::linalg::host_based::detail::extract_raw_pointer<float>(host_x), 1, 2,
66  viennacl::linalg::host_based::detail::extract_raw_pointer<float>(host_y), 0, 2);
67 
68  std::cout << " --- Host ---" << std::endl;
69  std::cout << "host_x: " << host_x << std::endl;
70  std::cout << "host_y: " << host_y << std::endl;
71 
78 #ifdef VIENNACL_WITH_CUDA
81 
82  ViennaCLCUDASswap(my_backend, half_size,
83  viennacl::cuda_arg(cuda_x), 0, 2,
84  viennacl::cuda_arg(cuda_y), 1, 2);
85 
86  std::cout << " --- CUDA ---" << std::endl;
87  std::cout << "cuda_x: " << cuda_x << std::endl;
88  std::cout << "cuda_y: " << cuda_y << std::endl;
89 #endif
90 
97 #ifdef VIENNACL_WITH_OPENCL
98  long context_id = 0;
101 
102  ViennaCLBackendSetOpenCLContextID(my_backend, static_cast<ViennaCLInt>(context_id));
103 
104  ViennaCLOpenCLSswap(my_backend, half_size,
105  viennacl::traits::opencl_handle(opencl_x).get(), 1, 2,
106  viennacl::traits::opencl_handle(opencl_y).get(), 1, 2);
107 
108  std::cout << " --- OpenCL ---" << std::endl;
109  std::cout << "opencl_x: " << opencl_x << std::endl;
110  std::cout << "opencl_y: " << opencl_y << std::endl;
111 #endif
112 
116  ViennaCLBackendDestroy(&my_backend);
117 
118  std::cout << "!!!! TUTORIAL COMPLETED SUCCESSFULLY !!!!" << std::endl;
119 
120  return EXIT_SUCCESS;
121 }
122 
VIENNACL_EXPORTED_FUNCTION ViennaCLStatus ViennaCLBackendCreate(ViennaCLBackend *backend)
Definition: backend.cpp:25
Generic backend for CUDA, OpenCL, host-based stuff.
int main()
Definition: bisect.cpp:91
VIENNACL_EXPORTED_FUNCTION ViennaCLStatus ViennaCLBackendSetOpenCLContextID(ViennaCLBackend backend, ViennaCLInt context_id)
Definition: backend.cpp:32
VIENNACL_EXPORTED_FUNCTION ViennaCLStatus ViennaCLHostSswap(ViennaCLBackend backend, ViennaCLInt n, float *x, ViennaCLInt offx, ViennaCLInt incx, float *y, ViennaCLInt offy, ViennaCLInt incy)
Definition: blas1_host.cpp:269
vcl_size_t size(VectorType const &vec)
Generic routine for obtaining the size of a vector (ViennaCL, uBLAS, etc.)
Definition: size.hpp:235
VIENNACL_EXPORTED_FUNCTION ViennaCLStatus ViennaCLBackendDestroy(ViennaCLBackend *backend)
Definition: backend.cpp:39
VIENNACL_EXPORTED_FUNCTION ViennaCLStatus ViennaCLCUDASswap(ViennaCLBackend backend, ViennaCLInt n, float *x, ViennaCLInt offx, ViennaCLInt incx, float *y, ViennaCLInt offy, ViennaCLInt incy)
int ViennaCLInt
Definition: viennacl.hpp:48
viennacl::context context(T const &t)
Returns an ID for the currently active memory domain of an object.
Definition: context.hpp:40
The vector type with operator-overloads and proxy classes is defined here. Linear algebra operations ...
Represents a vector consisting of scalars 's' only, i.e. v[i] = s for all i. To be used as an initial...
Definition: vector_def.hpp:87
NumericT * cuda_arg(scalar< NumericT > &obj)
Convenience helper function for extracting the CUDA handle from a ViennaCL scalar. Non-const version.
Definition: common.hpp:39
viennacl::ocl::context & get_context(long i)
Convenience function for returning the current context.
Definition: backend.hpp:225