ViennaCL - The Vienna Computing Library  1.6.0
Free open-source GPU-accelerated linear algebra and solver library.
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
compressed_matrix.hpp
Go to the documentation of this file.
1 #ifndef VIENNACL_COMPRESSED_MATRIX_HPP_
2 #define VIENNACL_COMPRESSED_MATRIX_HPP_
3 
4 /* =========================================================================
5  Copyright (c) 2010-2014, Institute for Microelectronics,
6  Institute for Analysis and Scientific Computing,
7  TU Wien.
8  Portions of this software are copyright by UChicago Argonne, LLC.
9 
10  -----------------
11  ViennaCL - The Vienna Computing Library
12  -----------------
13 
14  Project Head: Karl Rupp rupp@iue.tuwien.ac.at
15 
16  (A list of authors and contributors can be found in the PDF manual)
17 
18  License: MIT (X11), see file LICENSE in the base directory
19 ============================================================================= */
20 
25 #include <vector>
26 #include <list>
27 #include <map>
28 #include "viennacl/forwards.h"
29 #include "viennacl/vector.hpp"
30 
32 
33 #include "viennacl/tools/tools.hpp"
35 
36 #ifdef VIENNACL_WITH_UBLAS
37 #include <boost/numeric/ublas/matrix_sparse.hpp>
38 #endif
39 
40 namespace viennacl
41 {
42 namespace detail
43 {
44 
45  template<typename CPUMatrixT, typename NumericT, unsigned int AlignmentV>
46  void copy_impl(const CPUMatrixT & cpu_matrix,
48  vcl_size_t nonzeros)
49  {
50  assert( (gpu_matrix.size1() == 0 || viennacl::traits::size1(cpu_matrix) == gpu_matrix.size1()) && bool("Size mismatch") );
51  assert( (gpu_matrix.size2() == 0 || viennacl::traits::size2(cpu_matrix) == gpu_matrix.size2()) && bool("Size mismatch") );
52 
53  viennacl::backend::typesafe_host_array<unsigned int> row_buffer(gpu_matrix.handle1(), cpu_matrix.size1() + 1);
54  viennacl::backend::typesafe_host_array<unsigned int> col_buffer(gpu_matrix.handle2(), nonzeros);
55  std::vector<NumericT> elements(nonzeros);
56 
57  vcl_size_t row_index = 0;
58  vcl_size_t data_index = 0;
59 
60  for (typename CPUMatrixT::const_iterator1 row_it = cpu_matrix.begin1();
61  row_it != cpu_matrix.end1();
62  ++row_it)
63  {
64  row_buffer.set(row_index, data_index);
65  ++row_index;
66 
67  for (typename CPUMatrixT::const_iterator2 col_it = row_it.begin();
68  col_it != row_it.end();
69  ++col_it)
70  {
71  col_buffer.set(data_index, col_it.index2());
72  elements[data_index] = *col_it;
73  ++data_index;
74  }
75  data_index = viennacl::tools::align_to_multiple<vcl_size_t>(data_index, AlignmentV); //take care of alignment
76  }
77  row_buffer.set(row_index, data_index);
78 
79  gpu_matrix.set(row_buffer.get(),
80  col_buffer.get(),
81  &elements[0],
82  cpu_matrix.size1(),
83  cpu_matrix.size2(),
84  nonzeros);
85  }
86 }
87 
88 //provide copy-operation:
103 template<typename CPUMatrixT, typename NumericT, unsigned int AlignmentV>
104 void copy(const CPUMatrixT & cpu_matrix,
106 {
107  if ( cpu_matrix.size1() > 0 && cpu_matrix.size2() > 0 )
108  {
109  //determine nonzeros:
110  vcl_size_t num_entries = 0;
111  for (typename CPUMatrixT::const_iterator1 row_it = cpu_matrix.begin1();
112  row_it != cpu_matrix.end1();
113  ++row_it)
114  {
115  vcl_size_t entries_per_row = 0;
116  for (typename CPUMatrixT::const_iterator2 col_it = row_it.begin();
117  col_it != row_it.end();
118  ++col_it)
119  {
120  ++entries_per_row;
121  }
122  num_entries += viennacl::tools::align_to_multiple<vcl_size_t>(entries_per_row, AlignmentV);
123  }
124 
125  if (num_entries == 0) //we copy an empty matrix
126  num_entries = 1;
127 
128  //set up matrix entries:
129  viennacl::detail::copy_impl(cpu_matrix, gpu_matrix, num_entries);
130  }
131 }
132 
133 
134 //adapted for std::vector< std::map < > > argument:
140 template<typename SizeT, typename NumericT, unsigned int AlignmentV>
141 void copy(const std::vector< std::map<SizeT, NumericT> > & cpu_matrix,
143 {
144  vcl_size_t nonzeros = 0;
145  vcl_size_t max_col = 0;
146  for (vcl_size_t i=0; i<cpu_matrix.size(); ++i)
147  {
148  if (cpu_matrix[i].size() > 0)
149  nonzeros += ((cpu_matrix[i].size() - 1) / AlignmentV + 1) * AlignmentV;
150  if (cpu_matrix[i].size() > 0)
151  max_col = std::max<vcl_size_t>(max_col, (cpu_matrix[i].rbegin())->first);
152  }
153 
154  viennacl::detail::copy_impl(tools::const_sparse_matrix_adapter<NumericT, SizeT>(cpu_matrix, cpu_matrix.size(), max_col + 1),
155  gpu_matrix,
156  nonzeros);
157 }
158 
159 #ifdef VIENNACL_WITH_UBLAS
160 template<typename ScalarType, typename F, vcl_size_t IB, typename IA, typename TA>
161 void copy(const boost::numeric::ublas::compressed_matrix<ScalarType, F, IB, IA, TA> & ublas_matrix,
163 {
164  assert( (gpu_matrix.size1() == 0 || viennacl::traits::size1(ublas_matrix) == gpu_matrix.size1()) && bool("Size mismatch") );
165  assert( (gpu_matrix.size2() == 0 || viennacl::traits::size2(ublas_matrix) == gpu_matrix.size2()) && bool("Size mismatch") );
166 
167  //we just need to copy the CSR arrays:
168  viennacl::backend::typesafe_host_array<unsigned int> row_buffer(gpu_matrix.handle1(), ublas_matrix.size1() + 1);
169  for (vcl_size_t i=0; i<=ublas_matrix.size1(); ++i)
170  row_buffer.set(i, ublas_matrix.index1_data()[i]);
171 
172  viennacl::backend::typesafe_host_array<unsigned int> col_buffer(gpu_matrix.handle2(), ublas_matrix.nnz());
173  for (vcl_size_t i=0; i<ublas_matrix.nnz(); ++i)
174  col_buffer.set(i, ublas_matrix.index2_data()[i]);
175 
176  gpu_matrix.set(row_buffer.get(),
177  col_buffer.get(),
178  &(ublas_matrix.value_data()[0]),
179  ublas_matrix.size1(),
180  ublas_matrix.size2(),
181  ublas_matrix.nnz());
182 
183 }
184 #endif
185 
186 #ifdef VIENNACL_WITH_EIGEN
187 template<typename NumericT, int flags, unsigned int AlignmentV>
188 void copy(const Eigen::SparseMatrix<NumericT, flags> & eigen_matrix,
189  compressed_matrix<NumericT, AlignmentV> & gpu_matrix)
190 {
191  assert( (gpu_matrix.size1() == 0 || static_cast<vcl_size_t>(eigen_matrix.rows()) == gpu_matrix.size1()) && bool("Size mismatch") );
192  assert( (gpu_matrix.size2() == 0 || static_cast<vcl_size_t>(eigen_matrix.cols()) == gpu_matrix.size2()) && bool("Size mismatch") );
193 
194  std::vector< std::map<unsigned int, NumericT> > stl_matrix(eigen_matrix.rows());
195 
196  for (int k=0; k < eigen_matrix.outerSize(); ++k)
197  for (typename Eigen::SparseMatrix<NumericT, flags>::InnerIterator it(eigen_matrix, k); it; ++it)
198  stl_matrix[it.row()][it.col()] = it.value();
199 
200  copy(tools::const_sparse_matrix_adapter<NumericT>(stl_matrix, eigen_matrix.rows(), eigen_matrix.cols()), gpu_matrix);
201 }
202 #endif
203 
204 
205 #ifdef VIENNACL_WITH_MTL4
206 template<typename NumericT, unsigned int AlignmentV>
207 void copy(const mtl::compressed2D<NumericT> & cpu_matrix,
208  compressed_matrix<NumericT, AlignmentV> & gpu_matrix)
209 {
210  assert( (gpu_matrix.size1() == 0 || static_cast<vcl_size_t>(cpu_matrix.num_rows()) == gpu_matrix.size1()) && bool("Size mismatch") );
211  assert( (gpu_matrix.size2() == 0 || static_cast<vcl_size_t>(cpu_matrix.num_cols()) == gpu_matrix.size2()) && bool("Size mismatch") );
212 
213  typedef mtl::compressed2D<NumericT> MatrixType;
214 
215  std::vector< std::map<unsigned int, NumericT> > stl_matrix(cpu_matrix.num_rows());
216 
217  using mtl::traits::range_generator;
219 
220  // Choose between row and column traversal
221  typedef typename min<range_generator<mtl::tag::row, MatrixType>,
222  range_generator<mtl::tag::col, MatrixType> >::type range_type;
223  range_type my_range;
224 
225  // Type of outer cursor
226  typedef typename range_type::type c_type;
227  // Type of inner cursor
228  typedef typename mtl::traits::range_generator<mtl::tag::nz, c_type>::type ic_type;
229 
230  // Define the property maps
231  typename mtl::traits::row<MatrixType>::type row(cpu_matrix);
232  typename mtl::traits::col<MatrixType>::type col(cpu_matrix);
233  typename mtl::traits::const_value<MatrixType>::type value(cpu_matrix);
234 
235  // Now iterate over the matrix
236  for (c_type cursor(my_range.begin(cpu_matrix)), cend(my_range.end(cpu_matrix)); cursor != cend; ++cursor)
237  for (ic_type icursor(mtl::begin<mtl::tag::nz>(cursor)), icend(mtl::end<mtl::tag::nz>(cursor)); icursor != icend; ++icursor)
238  stl_matrix[row(*icursor)][col(*icursor)] = value(*icursor);
239 
240  copy(tools::const_sparse_matrix_adapter<NumericT>(stl_matrix, cpu_matrix.num_rows(), cpu_matrix.num_cols()), gpu_matrix);
241 }
242 #endif
243 
244 
245 
246 
247 
248 
249 
250 //
251 // gpu to cpu:
252 //
262 template<typename CPUMatrixT, typename NumericT, unsigned int AlignmentV>
264  CPUMatrixT & cpu_matrix )
265 {
266  assert( (viennacl::traits::size1(cpu_matrix) == gpu_matrix.size1()) && bool("Size mismatch") );
267  assert( (viennacl::traits::size2(cpu_matrix) == gpu_matrix.size2()) && bool("Size mismatch") );
268 
269  if ( gpu_matrix.size1() > 0 && gpu_matrix.size2() > 0 )
270  {
271  //get raw data from memory:
272  viennacl::backend::typesafe_host_array<unsigned int> row_buffer(gpu_matrix.handle1(), cpu_matrix.size1() + 1);
273  viennacl::backend::typesafe_host_array<unsigned int> col_buffer(gpu_matrix.handle2(), gpu_matrix.nnz());
274  std::vector<NumericT> elements(gpu_matrix.nnz());
275 
276  //std::cout << "GPU->CPU, nonzeros: " << gpu_matrix.nnz() << std::endl;
277 
278  viennacl::backend::memory_read(gpu_matrix.handle1(), 0, row_buffer.raw_size(), row_buffer.get());
279  viennacl::backend::memory_read(gpu_matrix.handle2(), 0, col_buffer.raw_size(), col_buffer.get());
280  viennacl::backend::memory_read(gpu_matrix.handle(), 0, sizeof(NumericT)* gpu_matrix.nnz(), &(elements[0]));
281 
282  //fill the cpu_matrix:
283  vcl_size_t data_index = 0;
284  for (vcl_size_t row = 1; row <= gpu_matrix.size1(); ++row)
285  {
286  while (data_index < row_buffer[row])
287  {
288  if (col_buffer[data_index] >= gpu_matrix.size2())
289  {
290  std::cerr << "ViennaCL encountered invalid data at colbuffer[" << data_index << "]: " << col_buffer[data_index] << std::endl;
291  return;
292  }
293 
294  if (std::fabs(elements[data_index]) > static_cast<NumericT>(0))
295  cpu_matrix(row-1, static_cast<vcl_size_t>(col_buffer[data_index])) = elements[data_index];
296  ++data_index;
297  }
298  }
299  }
300 }
301 
302 
308 template<typename NumericT, unsigned int AlignmentV>
310  std::vector< std::map<unsigned int, NumericT> > & cpu_matrix)
311 {
312  tools::sparse_matrix_adapter<NumericT> temp(cpu_matrix, cpu_matrix.size(), cpu_matrix.size());
313  copy(gpu_matrix, temp);
314 }
315 
316 #ifdef VIENNACL_WITH_UBLAS
317 template<typename ScalarType, unsigned int AlignmentV, typename F, vcl_size_t IB, typename IA, typename TA>
319  boost::numeric::ublas::compressed_matrix<ScalarType> & ublas_matrix)
320 {
321  assert( (viennacl::traits::size1(ublas_matrix) == gpu_matrix.size1()) && bool("Size mismatch") );
322  assert( (viennacl::traits::size2(ublas_matrix) == gpu_matrix.size2()) && bool("Size mismatch") );
323 
324  viennacl::backend::typesafe_host_array<unsigned int> row_buffer(gpu_matrix.handle1(), gpu_matrix.size1() + 1);
325  viennacl::backend::typesafe_host_array<unsigned int> col_buffer(gpu_matrix.handle2(), gpu_matrix.nnz());
326 
327  viennacl::backend::memory_read(gpu_matrix.handle1(), 0, row_buffer.raw_size(), row_buffer.get());
328  viennacl::backend::memory_read(gpu_matrix.handle2(), 0, col_buffer.raw_size(), col_buffer.get());
329 
330  ublas_matrix.clear();
331  ublas_matrix.reserve(gpu_matrix.nnz());
332 
333  ublas_matrix.set_filled(gpu_matrix.size1() + 1, gpu_matrix.nnz());
334 
335  for (vcl_size_t i=0; i<ublas_matrix.size1() + 1; ++i)
336  ublas_matrix.index1_data()[i] = row_buffer[i];
337 
338  for (vcl_size_t i=0; i<ublas_matrix.nnz(); ++i)
339  ublas_matrix.index2_data()[i] = col_buffer[i];
340 
341  viennacl::backend::memory_read(gpu_matrix.handle(), 0, sizeof(ScalarType) * gpu_matrix.nnz(), &(ublas_matrix.value_data()[0]));
342 
343 }
344 #endif
345 
346 #ifdef VIENNACL_WITH_EIGEN
347 template<typename NumericT, int flags, unsigned int AlignmentV>
348 void copy(compressed_matrix<NumericT, AlignmentV> & gpu_matrix,
349  Eigen::SparseMatrix<NumericT, flags> & eigen_matrix)
350 {
351  assert( (static_cast<vcl_size_t>(eigen_matrix.rows()) == gpu_matrix.size1()) && bool("Size mismatch") );
352  assert( (static_cast<vcl_size_t>(eigen_matrix.cols()) == gpu_matrix.size2()) && bool("Size mismatch") );
353 
354  if ( gpu_matrix.size1() > 0 && gpu_matrix.size2() > 0 )
355  {
356  //get raw data from memory:
357  viennacl::backend::typesafe_host_array<unsigned int> row_buffer(gpu_matrix.handle1(), gpu_matrix.size1() + 1);
358  viennacl::backend::typesafe_host_array<unsigned int> col_buffer(gpu_matrix.handle2(), gpu_matrix.nnz());
359  std::vector<NumericT> elements(gpu_matrix.nnz());
360 
361  viennacl::backend::memory_read(gpu_matrix.handle1(), 0, row_buffer.raw_size(), row_buffer.get());
362  viennacl::backend::memory_read(gpu_matrix.handle2(), 0, col_buffer.raw_size(), col_buffer.get());
363  viennacl::backend::memory_read(gpu_matrix.handle(), 0, sizeof(NumericT)* gpu_matrix.nnz(), &(elements[0]));
364 
365  eigen_matrix.setZero();
366  vcl_size_t data_index = 0;
367  for (vcl_size_t row = 1; row <= gpu_matrix.size1(); ++row)
368  {
369  while (data_index < row_buffer[row])
370  {
371  assert(col_buffer[data_index] < gpu_matrix.size2() && bool("ViennaCL encountered invalid data at col_buffer"));
372  if (elements[data_index] != static_cast<NumericT>(0.0))
373  eigen_matrix.insert(row-1, col_buffer[data_index]) = elements[data_index];
374  ++data_index;
375  }
376  }
377  }
378 }
379 #endif
380 
381 
382 
383 #ifdef VIENNACL_WITH_MTL4
384 template<typename NumericT, unsigned int AlignmentV>
385 void copy(compressed_matrix<NumericT, AlignmentV> & gpu_matrix,
386  mtl::compressed2D<NumericT> & mtl4_matrix)
387 {
388  assert( (static_cast<vcl_size_t>(mtl4_matrix.num_rows()) == gpu_matrix.size1()) && bool("Size mismatch") );
389  assert( (static_cast<vcl_size_t>(mtl4_matrix.num_cols()) == gpu_matrix.size2()) && bool("Size mismatch") );
390 
391  if ( gpu_matrix.size1() > 0 && gpu_matrix.size2() > 0 )
392  {
393 
394  //get raw data from memory:
395  viennacl::backend::typesafe_host_array<unsigned int> row_buffer(gpu_matrix.handle1(), gpu_matrix.size1() + 1);
396  viennacl::backend::typesafe_host_array<unsigned int> col_buffer(gpu_matrix.handle2(), gpu_matrix.nnz());
397  std::vector<NumericT> elements(gpu_matrix.nnz());
398 
399  viennacl::backend::memory_read(gpu_matrix.handle1(), 0, row_buffer.raw_size(), row_buffer.get());
400  viennacl::backend::memory_read(gpu_matrix.handle2(), 0, col_buffer.raw_size(), col_buffer.get());
401  viennacl::backend::memory_read(gpu_matrix.handle(), 0, sizeof(NumericT)* gpu_matrix.nnz(), &(elements[0]));
402 
403  //set_to_zero(mtl4_matrix);
404  //mtl4_matrix.change_dim(gpu_matrix.size1(), gpu_matrix.size2());
405 
406  mtl::matrix::inserter< mtl::compressed2D<NumericT> > ins(mtl4_matrix);
407  vcl_size_t data_index = 0;
408  for (vcl_size_t row = 1; row <= gpu_matrix.size1(); ++row)
409  {
410  while (data_index < row_buffer[row])
411  {
412  assert(col_buffer[data_index] < gpu_matrix.size2() && bool("ViennaCL encountered invalid data at col_buffer"));
413  if (elements[data_index] != static_cast<NumericT>(0.0))
414  ins(row-1, col_buffer[data_index]) << typename mtl::Collection< mtl::compressed2D<NumericT> >::value_type(elements[data_index]);
415  ++data_index;
416  }
417  }
418  }
419 }
420 #endif
421 
422 
423 
424 
425 
427 
432 template<class NumericT, unsigned int AlignmentV /* see VCLForwards.h */>
434 {
435 public:
439 
441  compressed_matrix() : rows_(0), cols_(0), nonzeros_(0) {}
442 
451  : rows_(rows), cols_(cols), nonzeros_(nonzeros)
452  {
453  row_buffer_.switch_active_handle_id(ctx.memory_type());
454  col_buffer_.switch_active_handle_id(ctx.memory_type());
455  elements_.switch_active_handle_id(ctx.memory_type());
456 
457 #ifdef VIENNACL_WITH_OPENCL
458  if (ctx.memory_type() == OPENCL_MEMORY)
459  {
460  row_buffer_.opencl_handle().context(ctx.opencl_context());
461  col_buffer_.opencl_handle().context(ctx.opencl_context());
462  elements_.opencl_handle().context(ctx.opencl_context());
463  }
464 #endif
465  if (rows > 0)
466  {
468  }
469  if (nonzeros > 0)
470  {
472  viennacl::backend::memory_create(elements_, sizeof(NumericT) * nonzeros, ctx);
473  }
474  }
475 
483  : rows_(rows), cols_(cols), nonzeros_(0)
484  {
485  row_buffer_.switch_active_handle_id(ctx.memory_type());
486  col_buffer_.switch_active_handle_id(ctx.memory_type());
487  elements_.switch_active_handle_id(ctx.memory_type());
488 
489 #ifdef VIENNACL_WITH_OPENCL
490  if (ctx.memory_type() == OPENCL_MEMORY)
491  {
492  row_buffer_.opencl_handle().context(ctx.opencl_context());
493  col_buffer_.opencl_handle().context(ctx.opencl_context());
494  elements_.opencl_handle().context(ctx.opencl_context());
495  }
496 #endif
497  if (rows > 0)
498  {
500  }
501  }
502 
503  explicit compressed_matrix(viennacl::context ctx) : rows_(0), cols_(0), nonzeros_(0)
504  {
505  row_buffer_.switch_active_handle_id(ctx.memory_type());
506  col_buffer_.switch_active_handle_id(ctx.memory_type());
507  elements_.switch_active_handle_id(ctx.memory_type());
508 
509 #ifdef VIENNACL_WITH_OPENCL
510  if (ctx.memory_type() == OPENCL_MEMORY)
511  {
512  row_buffer_.opencl_handle().context(ctx.opencl_context());
513  col_buffer_.opencl_handle().context(ctx.opencl_context());
514  elements_.opencl_handle().context(ctx.opencl_context());
515  }
516 #endif
517  }
518 
519 
520 #ifdef VIENNACL_WITH_OPENCL
521  explicit compressed_matrix(cl_mem mem_row_buffer, cl_mem mem_col_buffer, cl_mem mem_elements,
522  vcl_size_t rows, vcl_size_t cols, vcl_size_t nonzeros) :
523  rows_(rows), cols_(cols), nonzeros_(nonzeros)
524  {
526  row_buffer_.opencl_handle() = mem_row_buffer;
527  row_buffer_.opencl_handle().inc(); //prevents that the user-provided memory is deleted once the matrix object is destroyed.
528  row_buffer_.raw_size(sizeof(cl_uint) * (rows + 1));
529 
531  col_buffer_.opencl_handle() = mem_col_buffer;
532  col_buffer_.opencl_handle().inc(); //prevents that the user-provided memory is deleted once the matrix object is destroyed.
533  col_buffer_.raw_size(sizeof(cl_uint) * nonzeros);
534 
536  elements_.opencl_handle() = mem_elements;
537  elements_.opencl_handle().inc(); //prevents that the user-provided memory is deleted once the matrix object is destroyed.
538  elements_.raw_size(sizeof(NumericT) * nonzeros);
539  }
540 #endif
541 
542 
545  {
546  assert( (rows_ == 0 || rows_ == other.size1()) && bool("Size mismatch") );
547  assert( (cols_ == 0 || cols_ == other.size2()) && bool("Size mismatch") );
548 
549  rows_ = other.size1();
550  cols_ = other.size2();
551  nonzeros_ = other.nnz();
552 
553  viennacl::backend::typesafe_memory_copy<unsigned int>(other.row_buffer_, row_buffer_);
554  viennacl::backend::typesafe_memory_copy<unsigned int>(other.col_buffer_, col_buffer_);
555  viennacl::backend::typesafe_memory_copy<NumericT>(other.elements_, elements_);
556 
557  return *this;
558  }
559 
560 
570  void set(const void * row_jumper,
571  const void * col_buffer,
572  const NumericT * elements,
573  vcl_size_t rows,
574  vcl_size_t cols,
575  vcl_size_t nonzeros)
576  {
577  assert( (rows > 0) && bool("Error in compressed_matrix::set(): Number of rows must be larger than zero!"));
578  assert( (cols > 0) && bool("Error in compressed_matrix::set(): Number of columns must be larger than zero!"));
579  assert( (nonzeros > 0) && bool("Error in compressed_matrix::set(): Number of nonzeros must be larger than zero!"));
580  //std::cout << "Setting memory: " << cols + 1 << ", " << nonzeros << std::endl;
581 
582  //row_buffer_.switch_active_handle_id(viennacl::backend::OPENCL_MEMORY);
584 
585  //col_buffer_.switch_active_handle_id(viennacl::backend::OPENCL_MEMORY);
587 
588  //elements_.switch_active_handle_id(viennacl::backend::OPENCL_MEMORY);
589  viennacl::backend::memory_create(elements_, sizeof(NumericT) * nonzeros, viennacl::traits::context(elements_), elements);
590 
591  nonzeros_ = nonzeros;
592  rows_ = rows;
593  cols_ = cols;
594  }
595 
597  void reserve(vcl_size_t new_nonzeros)
598  {
599  if (new_nonzeros > nonzeros_)
600  {
601  handle_type col_buffer_old;
602  handle_type elements_old;
603  viennacl::backend::memory_shallow_copy(col_buffer_, col_buffer_old);
604  viennacl::backend::memory_shallow_copy(elements_, elements_old);
605 
607  viennacl::backend::memory_create(col_buffer_, size_deducer.element_size() * new_nonzeros, viennacl::traits::context(col_buffer_));
608  viennacl::backend::memory_create(elements_, sizeof(NumericT) * new_nonzeros, viennacl::traits::context(elements_));
609 
610  viennacl::backend::memory_copy(col_buffer_old, col_buffer_, 0, 0, size_deducer.element_size() * nonzeros_);
611  viennacl::backend::memory_copy(elements_old, elements_, 0, 0, sizeof(NumericT)* nonzeros_);
612 
613  nonzeros_ = new_nonzeros;
614  }
615  }
616 
623  void resize(vcl_size_t new_size1, vcl_size_t new_size2, bool preserve = true)
624  {
625  assert(new_size1 > 0 && new_size2 > 0 && bool("Cannot resize to zero size!"));
626 
627  if (new_size1 != rows_ || new_size2 != cols_)
628  {
629  std::vector<std::map<unsigned int, NumericT> > stl_sparse_matrix;
630  if (rows_ > 0)
631  {
632  if (preserve)
633  {
634  stl_sparse_matrix.resize(rows_);
635  viennacl::copy(*this, stl_sparse_matrix);
636  } else
637  stl_sparse_matrix[0][0] = 0;
638  } else {
639  stl_sparse_matrix.resize(new_size1);
640  stl_sparse_matrix[0][0] = 0; //enforces nonzero array sizes if matrix was initially empty
641  }
642 
643  stl_sparse_matrix.resize(new_size1);
644 
645  //discard entries with column index larger than new_size2
646  if (new_size2 < cols_ && rows_ > 0)
647  {
648  for (vcl_size_t i=0; i<stl_sparse_matrix.size(); ++i)
649  {
650  std::list<unsigned int> to_delete;
651  for (typename std::map<unsigned int, NumericT>::iterator it = stl_sparse_matrix[i].begin();
652  it != stl_sparse_matrix[i].end();
653  ++it)
654  {
655  if (it->first >= new_size2)
656  to_delete.push_back(it->first);
657  }
658 
659  for (std::list<unsigned int>::iterator it = to_delete.begin(); it != to_delete.end(); ++it)
660  stl_sparse_matrix[i].erase(*it);
661  }
662  }
663 
664  viennacl::copy(stl_sparse_matrix, *this);
665 
666  rows_ = new_size1;
667  cols_ = new_size2;
668  }
669  }
670 
672  void clear()
673  {
674  viennacl::backend::typesafe_host_array<unsigned int> host_row_buffer(row_buffer_, rows_ + 1);
675  viennacl::backend::typesafe_host_array<unsigned int> host_col_buffer(col_buffer_, 1);
676  std::vector<NumericT> host_elements(1);
677 
678  viennacl::backend::memory_create(row_buffer_, host_row_buffer.element_size() * (rows_ + 1), viennacl::traits::context(row_buffer_), host_row_buffer.get());
679  viennacl::backend::memory_create(col_buffer_, host_col_buffer.element_size() * 1, viennacl::traits::context(col_buffer_), host_col_buffer.get());
680  viennacl::backend::memory_create(elements_, sizeof(NumericT) * 1, viennacl::traits::context(elements_), &(host_elements[0]));
681 
682  nonzeros_ = 0;
683  }
684 
687  {
688  assert( (i < rows_) && (j < cols_) && bool("compressed_matrix access out of bounds!"));
689 
690  vcl_size_t index = element_index(i, j);
691 
692  // check for element in sparsity pattern
693  if (index < nonzeros_)
694  return entry_proxy<NumericT>(index, elements_);
695 
696  // Element not found. Copying required. Very slow, but direct entry manipulation is painful anyway...
697  std::vector< std::map<unsigned int, NumericT> > cpu_backup(rows_);
698  tools::sparse_matrix_adapter<NumericT> adapted_cpu_backup(cpu_backup, rows_, cols_);
699  viennacl::copy(*this, adapted_cpu_backup);
700  cpu_backup[i][static_cast<unsigned int>(j)] = 0.0;
701  viennacl::copy(adapted_cpu_backup, *this);
702 
703  index = element_index(i, j);
704 
705  assert(index < nonzeros_);
706 
707  return entry_proxy<NumericT>(index, elements_);
708  }
709 
711  const vcl_size_t & size1() const { return rows_; }
713  const vcl_size_t & size2() const { return cols_; }
715  const vcl_size_t & nnz() const { return nonzeros_; }
716 
718  const handle_type & handle1() const { return row_buffer_; }
720  const handle_type & handle2() const { return col_buffer_; }
722  const handle_type & handle() const { return elements_; }
723 
725  handle_type & handle1() { return row_buffer_; }
727  handle_type & handle2() { return col_buffer_; }
729  handle_type & handle() { return elements_; }
730 
732  {
733  viennacl::backend::switch_memory_context<unsigned int>(row_buffer_, new_ctx);
734  viennacl::backend::switch_memory_context<unsigned int>(col_buffer_, new_ctx);
735  viennacl::backend::switch_memory_context<NumericT>(elements_, new_ctx);
736  }
737 
739  {
740  return row_buffer_.get_active_handle_id();
741  }
742 
743 private:
744 
745  vcl_size_t element_index(vcl_size_t i, vcl_size_t j)
746  {
747  //read row indices
748  viennacl::backend::typesafe_host_array<unsigned int> row_indices(row_buffer_, 2);
749  viennacl::backend::memory_read(row_buffer_, row_indices.element_size()*i, row_indices.element_size()*2, row_indices.get());
750 
751  //get column indices for row i:
752  viennacl::backend::typesafe_host_array<unsigned int> col_indices(col_buffer_, row_indices[1] - row_indices[0]);
753  viennacl::backend::memory_read(col_buffer_, col_indices.element_size()*row_indices[0], row_indices.element_size()*col_indices.size(), col_indices.get());
754 
755  for (vcl_size_t k=0; k<col_indices.size(); ++k)
756  {
757  if (col_indices[k] == j)
758  return row_indices[0] + k;
759  }
760 
761  // if not found, return index past the end of the matrix (cf. matrix.end() in the spirit of the STL)
762  return nonzeros_;
763  }
764 
765  // /** @brief Copy constructor is by now not available. */
766  //compressed_matrix(compressed_matrix const &);
767 
768 
769  vcl_size_t rows_;
770  vcl_size_t cols_;
771  vcl_size_t nonzeros_;
772  handle_type row_buffer_;
773  handle_type col_buffer_;
774  handle_type elements_;
775 };
776 
777 
778 
779 //
780 // Specify available operations:
781 //
782 
785 namespace linalg
786 {
787 namespace detail
788 {
789  // x = A * y
790  template<typename T, unsigned int A>
791  struct op_executor<vector_base<T>, op_assign, vector_expression<const compressed_matrix<T, A>, const vector_base<T>, op_prod> >
792  {
793  static void apply(vector_base<T> & lhs, vector_expression<const compressed_matrix<T, A>, const vector_base<T>, op_prod> const & rhs)
794  {
795  // check for the special case x = A * x
796  if (viennacl::traits::handle(lhs) == viennacl::traits::handle(rhs.rhs()))
797  {
798  viennacl::vector<T> temp(lhs);
799  viennacl::linalg::prod_impl(rhs.lhs(), rhs.rhs(), temp);
800  lhs = temp;
801  }
802  else
803  viennacl::linalg::prod_impl(rhs.lhs(), rhs.rhs(), lhs);
804  }
805  };
806 
807  template<typename T, unsigned int A>
808  struct op_executor<vector_base<T>, op_inplace_add, vector_expression<const compressed_matrix<T, A>, const vector_base<T>, op_prod> >
809  {
810  static void apply(vector_base<T> & lhs, vector_expression<const compressed_matrix<T, A>, const vector_base<T>, op_prod> const & rhs)
811  {
812  viennacl::vector<T> temp(lhs);
813  viennacl::linalg::prod_impl(rhs.lhs(), rhs.rhs(), temp);
814  lhs += temp;
815  }
816  };
817 
818  template<typename T, unsigned int A>
819  struct op_executor<vector_base<T>, op_inplace_sub, vector_expression<const compressed_matrix<T, A>, const vector_base<T>, op_prod> >
820  {
821  static void apply(vector_base<T> & lhs, vector_expression<const compressed_matrix<T, A>, const vector_base<T>, op_prod> const & rhs)
822  {
823  viennacl::vector<T> temp(lhs);
824  viennacl::linalg::prod_impl(rhs.lhs(), rhs.rhs(), temp);
825  lhs -= temp;
826  }
827  };
828 
829 
830  // x = A * vec_op
831  template<typename T, unsigned int A, typename LHS, typename RHS, typename OP>
832  struct op_executor<vector_base<T>, op_assign, vector_expression<const compressed_matrix<T, A>, const vector_expression<const LHS, const RHS, OP>, op_prod> >
833  {
834  static void apply(vector_base<T> & lhs, vector_expression<const compressed_matrix<T, A>, const vector_expression<const LHS, const RHS, OP>, op_prod> const & rhs)
835  {
836  viennacl::vector<T> temp(rhs.rhs(), viennacl::traits::context(rhs));
837  viennacl::linalg::prod_impl(rhs.lhs(), temp, lhs);
838  }
839  };
840 
841  // x = A * vec_op
842  template<typename T, unsigned int A, typename LHS, typename RHS, typename OP>
843  struct op_executor<vector_base<T>, op_inplace_add, vector_expression<const compressed_matrix<T, A>, vector_expression<const LHS, const RHS, OP>, op_prod> >
844  {
845  static void apply(vector_base<T> & lhs, vector_expression<const compressed_matrix<T, A>, vector_expression<const LHS, const RHS, OP>, op_prod> const & rhs)
846  {
847  viennacl::vector<T> temp(rhs.rhs(), viennacl::traits::context(rhs));
848  viennacl::vector<T> temp_result(lhs);
849  viennacl::linalg::prod_impl(rhs.lhs(), temp, temp_result);
850  lhs += temp_result;
851  }
852  };
853 
854  // x = A * vec_op
855  template<typename T, unsigned int A, typename LHS, typename RHS, typename OP>
856  struct op_executor<vector_base<T>, op_inplace_sub, vector_expression<const compressed_matrix<T, A>, const vector_expression<const LHS, const RHS, OP>, op_prod> >
857  {
858  static void apply(vector_base<T> & lhs, vector_expression<const compressed_matrix<T, A>, const vector_expression<const LHS, const RHS, OP>, op_prod> const & rhs)
859  {
860  viennacl::vector<T> temp(rhs.rhs(), viennacl::traits::context(rhs));
861  viennacl::vector<T> temp_result(lhs);
862  viennacl::linalg::prod_impl(rhs.lhs(), temp, temp_result);
863  lhs -= temp_result;
864  }
865  };
866 
867 } // namespace detail
868 } // namespace linalg
870 }
871 
872 #endif
const vcl_size_t & size2() const
Returns the number of columns.
Helper class implementing an array on the host. Default case: No conversion necessary.
Definition: util.hpp:92
vcl_size_t element_size() const
Definition: util.hpp:112
This class represents a single scalar value on the GPU and behaves mostly like a built-in scalar type...
Definition: forwards.h:226
const vcl_size_t & size1() const
Returns the number of rows.
Various little tools used here and there in ViennaCL.
void switch_memory_context(viennacl::context new_ctx)
vcl_size_t size1(MatrixType const &mat)
Generic routine for obtaining the number of rows of a matrix (ViennaCL, uBLAS, etc.)
Definition: size.hpp:216
A proxy class for entries in a vector.
This file provides the forward declarations for the main types used within ViennaCL.
void memory_read(mem_handle const &src_buffer, vcl_size_t src_offset, vcl_size_t bytes_to_read, void *ptr, bool async=false)
Reads data from a buffer back to main RAM.
Definition: memory.hpp:261
result_of::size_type< MatrixType >::type size2(MatrixType const &mat)
Generic routine for obtaining the number of columns of a matrix (ViennaCL, uBLAS, etc...
Definition: size.hpp:245
const handle_type & handle() const
Returns the OpenCL handle to the matrix entry array.
const handle_type & handle1() const
Returns the OpenCL handle to the row index array.
const vcl_size_t & nnz() const
Returns the number of nonzero entries.
vcl_size_t element_size(memory_types)
Definition: memory.hpp:299
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:144
entry_proxy< NumericT > operator()(vcl_size_t i, vcl_size_t j)
Returns a reference to the (i,j)-th entry of the sparse matrix. If (i,j) does not exist (zero)...
handle_type & handle()
Returns the OpenCL handle to the matrix entry array.
Implementations of operations using sparse matrices.
const handle_type & handle2() const
Returns the OpenCL handle to the column index array.
void clear()
Resets all entries in the matrix back to zero without changing the matrix size. Resets the sparsity p...
void reserve(vcl_size_t new_nonzeros)
Allocate memory for the supplied number of nonzeros in the matrix. Old values are preserved...
Adapts a constant sparse matrix type made up from std::vector > to basic ub...
Definition: adapter.hpp:183
scalar< typename viennacl::tools::CHECK_SCALAR_TEMPLATE_ARGUMENT< NumericT >::ResultType > value_type
std::size_t vcl_size_t
Definition: forwards.h:74
compressed_matrix(vcl_size_t rows, vcl_size_t cols, viennacl::context ctx)
Construction of a compressed matrix with the supplied number of rows and columns. If the number of no...
vector_expression< const matrix_base< NumericT, F >, const unsigned int, op_row > row(const matrix_base< NumericT, F > &A, unsigned int i)
Definition: matrix.hpp:853
viennacl::memory_types memory_type() const
Definition: context.hpp:76
void copy_impl(const CPUMatrixT &cpu_matrix, compressed_compressed_matrix< NumericT > &gpu_matrix, vcl_size_t nonzero_rows, vcl_size_t nonzeros)
viennacl::memory_types memory_context() const
handle_type & handle1()
Returns the OpenCL handle to the row index array.
compressed_matrix(vcl_size_t rows, vcl_size_t cols, vcl_size_t nonzeros=0, viennacl::context ctx=viennacl::context())
Construction of a compressed matrix with the supplied number of rows and columns. If the number of no...
void switch_active_handle_id(memory_types new_id)
Switches the currently active handle. If no support for that backend is provided, an exception is thr...
Definition: mem_handle.hpp:121
void memory_copy(mem_handle const &src_buffer, mem_handle &dst_buffer, vcl_size_t src_offset, vcl_size_t dst_offset, vcl_size_t bytes_to_copy)
Copies 'bytes_to_copy' bytes from address 'src_buffer + src_offset' to memory starting at address 'ds...
Definition: memory.hpp:140
compressed_matrix(viennacl::context ctx)
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 ...
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) ...
void set(const void *row_jumper, const void *col_buffer, const NumericT *elements, vcl_size_t rows, vcl_size_t cols, vcl_size_t nonzeros)
Sets the row, column and value arrays of the compressed matrix.
void set(vcl_size_t index, U value)
Definition: util.hpp:115
viennacl::backend::mem_handle handle_type
float ScalarType
Definition: fft_1d.cpp:42
Main abstraction class for multiple memory domains. Represents a buffer in either main RAM...
Definition: mem_handle.hpp:89
vcl_size_t raw_size() const
Returns the number of bytes of the currently active buffer.
Definition: mem_handle.hpp:230
Adapts a non-const sparse matrix type made up from std::vector > to basic u...
Definition: adapter.hpp:357
A sparse square matrix in compressed sparse rows format.
void memory_create(mem_handle &handle, vcl_size_t size_in_bytes, viennacl::context const &ctx, const void *host_ptr=NULL)
Creates an array of the specified size. If the second argument is provided, the buffer is initialized...
Definition: memory.hpp:87
T min(const T &lhs, const T &rhs)
Minimum.
Definition: util.hpp:45
void prod_impl(const matrix_base< NumericT > &mat, const vector_base< NumericT > &vec, vector_base< NumericT > &result)
Carries out matrix-vector multiplication.
viennacl::backend::mem_handle & handle(T &obj)
Returns the generic memory handle of an object. Non-const version.
Definition: handle.hpp:41
memory_types
Definition: forwards.h:344
A proxy class for a single element of a vector or matrix. This proxy should not be noticed by end-use...
Definition: forwards.h:232
compressed_matrix()
Default construction of a compressed matrix. No memory is allocated.
void resize(vcl_size_t new_size1, vcl_size_t new_size2, bool preserve=true)
Resize the matrix.
void memory_shallow_copy(mem_handle const &src_buffer, mem_handle &dst_buffer)
A 'shallow' copy operation from an initialized buffer to an uninitialized buffer. The uninitialized b...
Definition: memory.hpp:177
compressed_matrix & operator=(compressed_matrix const &other)
Assignment a compressed matrix from possibly another memory domain.
handle_type & handle2()
Returns the OpenCL handle to the column index array.
memory_types get_active_handle_id() const
Returns an ID for the currently active memory buffer. Other memory buffers might contain old or no da...
Definition: mem_handle.hpp:118