ViennaCL - The Vienna Computing Library  1.7.0
Free open-source GPU-accelerated linear algebra and solver library.
multithreaded_cg.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 
25 #ifndef VIENNACL_WITH_OPENCL
26  #define VIENNACL_WITH_OPENCL
27 #endif
28 
29 // include necessary system headers
30 #include <iostream>
31 
32 //
33 // ublas includes
34 //
35 #include <boost/numeric/ublas/io.hpp>
36 #include <boost/numeric/ublas/triangular.hpp>
37 #include <boost/numeric/ublas/matrix_sparse.hpp>
38 #include <boost/numeric/ublas/matrix.hpp>
39 #include <boost/numeric/ublas/matrix_proxy.hpp>
40 #include <boost/numeric/ublas/operation.hpp>
41 #include <boost/numeric/ublas/operation_sparse.hpp>
42 #include <boost/numeric/ublas/io.hpp>
43 #include <boost/numeric/ublas/lu.hpp>
44 
45 // Must be set if you want to use ViennaCL algorithms on ublas objects
46 #define VIENNACL_WITH_UBLAS 1
47 
48 
49 //include basic scalar and vector types of ViennaCL
50 #include "viennacl/scalar.hpp"
51 #include "viennacl/vector.hpp"
54 
55 #include "viennacl/ocl/device.hpp"
57 #include "viennacl/ocl/backend.hpp"
58 
59 //include the generic inner product functions of ViennaCL
61 #include "viennacl/linalg/cg.hpp"
62 
63 // Some helper functions for this tutorial:
64 #include "vector-io.hpp"
65 
66 using namespace boost::numeric;
67 
71 #include <boost/thread.hpp>
72 
77 template<typename NumericT>
78 class worker
79 {
80 public:
81  worker(std::size_t tid) : thread_id_(tid) {}
82 
86  void operator()()
87  {
91  ublas::vector<NumericT> rhs;
92  ublas::vector<NumericT> ref_result;
93  ublas::compressed_matrix<NumericT> ublas_matrix;
94 
98  if (!viennacl::io::read_matrix_market_file(ublas_matrix, "../examples/testdata/mat65k.mtx"))
99  {
100  std::cout << "Error reading Matrix file" << std::endl;
101  return;
102  }
103 
104  if (!readVectorFromFile("../examples/testdata/rhs65025.txt", rhs))
105  {
106  std::cout << "Error reading RHS file" << std::endl;
107  return;
108  }
109 
110  if (!readVectorFromFile("../examples/testdata/result65025.txt", ref_result))
111  {
112  std::cout << "Error reading Result file" << std::endl;
113  return;
114  }
115 
120  viennacl::context ctx(viennacl::ocl::get_context(static_cast<long>(thread_id_)));
121 
122  std::size_t vcl_size = rhs.size();
123  viennacl::compressed_matrix<NumericT> vcl_compressed_matrix(ctx);
124  viennacl::vector<NumericT> vcl_rhs(vcl_size, ctx);
125  viennacl::vector<NumericT> vcl_ref_result(vcl_size, ctx);
126 
127  viennacl::copy(rhs.begin(), rhs.end(), vcl_rhs.begin());
128  viennacl::copy(ref_result.begin(), ref_result.end(), vcl_ref_result.begin());
129 
130 
134  viennacl::copy(ublas_matrix, vcl_compressed_matrix);
135 
136  viennacl::vector<NumericT> vcl_result = viennacl::linalg::solve(vcl_compressed_matrix, vcl_rhs, viennacl::linalg::cg_tag());
137 
138  std::stringstream ss;
139  ss << "Result of thread " << thread_id_ << " on device " << viennacl::ocl::get_context(static_cast<long>(thread_id_)).devices()[0].name() << ": " << vcl_result[0] << ", should: " << ref_result[0] << std::endl;
140  message_ = ss.str();
141  }
142 
143  std::string message() const { return message_; }
144 
145 private:
146  std::string message_;
147  std::size_t thread_id_;
148 };
149 
153 int main()
154 {
155  //Change this type definition to double if your gpu supports that
156  typedef float ScalarType;
157 
158  if (viennacl::ocl::get_platforms().size() == 0)
159  {
160  std::cerr << "Error: No platform found!" << std::endl;
161  return EXIT_FAILURE;
162  }
163 
168  std::vector<viennacl::ocl::device> const & devices = pf.devices();
169 
170  // Set first device to first context:
171  viennacl::ocl::setup_context(0, devices[0]);
172 
173  // Set second device for second context (use the same device for the second context if only one device available):
174  if (devices.size() > 1)
175  viennacl::ocl::setup_context(1, devices[1]);
176  else
177  viennacl::ocl::setup_context(1, devices[0]);
178 
183  worker<ScalarType> work_functor0(0);
184  worker<ScalarType> work_functor1(1);
185  boost::thread worker_thread_0(boost::ref(work_functor0));
186  boost::thread worker_thread_1(boost::ref(work_functor1));
187 
188  worker_thread_0.join();
189  worker_thread_1.join();
190 
191  std::cout << work_functor0.message() << std::endl;
192  std::cout << work_functor1.message() << std::endl;
193 
197  std::cout << "!!!! TUTORIAL COMPLETED SUCCESSFULLY !!!!" << std::endl;
198 
199  return EXIT_SUCCESS;
200 }
201 
A reader and writer for the matrix market format is implemented here.
std::vector< platform > get_platforms()
Definition: platform.hpp:124
Generic interface for the l^2-norm. See viennacl/linalg/vector_operations.hpp for implementations...
Represents an OpenCL device within ViennaCL.
Implements a OpenCL platform within ViennaCL.
Wrapper class for an OpenCL platform.
Definition: platform.hpp:45
int main()
Definition: bisect.cpp:91
std::vector< device > devices(cl_device_type dtype=CL_DEVICE_TYPE_DEFAULT)
Returns the available devices of the supplied device type.
Definition: platform.hpp:91
VectorT solve(MatrixT const &matrix, VectorT const &rhs, bicgstab_tag const &tag, PreconditionerT const &precond)
Definition: bicgstab.hpp:496
Represents a generic 'context' similar to an OpenCL context, but is backend-agnostic and thus also su...
Definition: context.hpp:39
vcl_size_t size(VectorType const &vec)
Generic routine for obtaining the size of a vector (ViennaCL, uBLAS, etc.)
Definition: size.hpp:235
Implementation of the compressed_matrix class.
bool readVectorFromFile(const std::string &filename, VectorType &vec)
Definition: vector-io.hpp:104
The conjugate gradient method is implemented here.
Implementations of the OpenCL backend, where all contexts are stored in.
The vector type with operator-overloads and proxy classes is defined here. Linear algebra operations ...
void copy(std::vector< NumericT > &cpu_vec, circulant_matrix< NumericT, AlignmentV > &gpu_mat)
Copies a circulant matrix from the std::vector to the OpenCL device (either GPU or multi-core CPU) ...
float ScalarType
Definition: fft_1d.cpp:42
A tag for the conjugate gradient Used for supplying solver parameters and for dispatching the solve()...
Definition: cg.hpp:48
long read_matrix_market_file(MatrixT &mat, const char *file, long index_base=1)
Reads a sparse matrix from a file (MatrixMarket format)
viennacl::ocl::context & get_context(long i)
Convenience function for returning the current context.
Definition: backend.hpp:225
std::vector< viennacl::ocl::device > const & devices() const
Returns a vector with all devices in this context.
Definition: context.hpp:106
Implementation of the ViennaCL scalar class.
void setup_context(long i, std::vector< cl_device_id > const &devices)
Convenience function for setting devices for a context.
Definition: backend.hpp:231