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
size.hpp
Go to the documentation of this file.
1 #ifndef VIENNACL_TRAITS_SIZE_HPP_
2 #define VIENNACL_TRAITS_SIZE_HPP_
3 
4 /* =========================================================================
5  Copyright (c) 2010-2016, 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 manual)
17 
18  License: MIT (X11), see file LICENSE in the base directory
19 ============================================================================= */
20 
25 #include <string>
26 #include <fstream>
27 #include <sstream>
28 #include "viennacl/forwards.h"
31 
32 #ifdef VIENNACL_WITH_UBLAS
33 #include <boost/numeric/ublas/matrix_sparse.hpp>
34 #include <boost/numeric/ublas/matrix.hpp>
35 #endif
36 
37 #ifdef VIENNACL_WITH_ARMADILLO
38 #include <armadillo>
39 #endif
40 
41 #ifdef VIENNACL_WITH_EIGEN
42 #include <Eigen/Core>
43 #include <Eigen/Sparse>
44 #endif
45 
46 #ifdef VIENNACL_WITH_MTL4
47 #include <boost/numeric/mtl/mtl.hpp>
48 #endif
49 
50 #include <vector>
51 #include <map>
52 
53 namespace viennacl
54 {
55 namespace traits
56 {
57 
58 //
59 // Resize: Change the size of vectors and matrices
60 //
62 template<typename MatrixType>
63 void resize(MatrixType & matrix, vcl_size_t rows, vcl_size_t cols)
64 {
65  matrix.resize(rows, cols);
66 }
67 
69 template<typename VectorType>
70 void resize(VectorType & vec, vcl_size_t new_size)
71 {
72  vec.resize(new_size);
73 }
74 
76 #ifdef VIENNACL_WITH_UBLAS
77 //ublas needs separate treatment:
78 template<typename ScalarType>
79 void resize(boost::numeric::ublas::compressed_matrix<ScalarType> & matrix,
80  vcl_size_t rows,
81  vcl_size_t cols)
82 {
83  matrix.resize(rows, cols, false); //Note: omitting third parameter leads to compile time error (not implemented in ublas <= 1.42)
84 }
85 #endif
86 
87 
88 #ifdef VIENNACL_WITH_MTL4
89 template<typename ScalarType>
90 void resize(mtl::compressed2D<ScalarType> & matrix,
91  vcl_size_t rows,
92  vcl_size_t cols)
93 {
94  matrix.change_dim(rows, cols);
95 }
96 
97 template<typename ScalarType>
98 void resize(mtl::dense_vector<ScalarType> & vec,
99  vcl_size_t new_size)
100 {
101  vec.change_dim(new_size);
102 }
103 #endif
104 
105 #ifdef VIENNACL_WITH_ARMADILLO
106 template<typename NumericT>
107 inline void resize(arma::Mat<NumericT> & A,
108  vcl_size_t new_rows,
109  vcl_size_t new_cols)
110 {
111  A.resize(new_rows, new_cols);
112 }
113 
114 template<typename NumericT>
115 inline void resize(arma::SpMat<NumericT> & A,
116  vcl_size_t new_rows,
117  vcl_size_t new_cols)
118 {
119  A.set_size(new_rows, new_cols);
120 }
121 #endif
122 
123 #ifdef VIENNACL_WITH_EIGEN
124 template<typename NumericT, int Options>
125 inline void resize(Eigen::Matrix<NumericT, Eigen::Dynamic, Eigen::Dynamic, Options> & m,
126  vcl_size_t new_rows,
127  vcl_size_t new_cols)
128 {
129  m.resize(new_rows, new_cols);
130 }
131 
132 template<typename T, int options>
133 inline void resize(Eigen::SparseMatrix<T, options> & m,
134  vcl_size_t new_rows,
135  vcl_size_t new_cols)
136 {
137  m.resize(new_rows, new_cols);
138 }
139 
140 inline void resize(Eigen::VectorXf & v,
141  vcl_size_t new_size)
142 {
143  v.resize(new_size);
144 }
145 
146 inline void resize(Eigen::VectorXd & v,
147  vcl_size_t new_size)
148 {
149  v.resize(new_size);
150 }
151 #endif
152 
157 //
158 // size1: No. of rows for matrices
159 //
161 template<typename MatrixType>
163 size1(MatrixType const & mat) { return mat.size1(); }
164 
166 template<typename RowType>
168 size1(std::vector< RowType > const & mat) { return mat.size(); }
169 
170 #ifdef VIENNACL_WITH_ARMADILLO
171 template<typename NumericT>
172 inline vcl_size_t size1(arma::Mat<NumericT> const & A) { return A.n_rows; }
173 template<typename NumericT>
174 inline vcl_size_t size1(arma::SpMat<NumericT> const & A) { return A.n_rows; }
175 #endif
176 
177 #ifdef VIENNACL_WITH_EIGEN
178 template<typename NumericT, int Options>
179 vcl_size_t size1(Eigen::Matrix<NumericT, Eigen::Dynamic, Eigen::Dynamic, Options> const & m) { return static_cast<vcl_size_t>(m.rows()); }
180 template<typename NumericT, int Options>
181 vcl_size_t size1(Eigen::Map< Eigen::Matrix<NumericT, Eigen::Dynamic, Eigen::Dynamic, Options> > const & m) { return static_cast<vcl_size_t>(m.rows()); }
182 template<typename T, int options>
183 inline vcl_size_t size1(Eigen::SparseMatrix<T, options> & m) { return static_cast<vcl_size_t>(m.rows()); }
184 #endif
185 
186 #ifdef VIENNACL_WITH_MTL4
187 template<typename NumericT, typename T>
188 vcl_size_t size1(mtl::dense2D<NumericT, T> const & m) { return static_cast<vcl_size_t>(m.num_rows()); }
189 template<typename NumericT>
190 vcl_size_t size1(mtl::compressed2D<NumericT> const & m) { return static_cast<vcl_size_t>(m.num_rows()); }
191 #endif
192 
195 //
196 // size2: No. of columns for matrices
197 //
199 template<typename MatrixType>
200 typename result_of::size_type<MatrixType>::type
201 size2(MatrixType const & mat) { return mat.size2(); }
202 
204 template<typename RowType>
206 size2(std::vector< RowType > const & mat) { return mat[0].size(); }
207 
208 #ifdef VIENNACL_WITH_ARMADILLO
209 template<typename NumericT>
210 inline vcl_size_t size2(arma::Mat<NumericT> const & A) { return A.n_cols; }
211 template<typename NumericT>
212 inline vcl_size_t size2(arma::SpMat<NumericT> const & A) { return A.n_cols; }
213 #endif
214 
215 #ifdef VIENNACL_WITH_EIGEN
216 template<typename NumericT, int Options>
217 inline vcl_size_t size2(Eigen::Matrix<NumericT, Eigen::Dynamic, Eigen::Dynamic, Options> const & m) { return m.cols(); }
218 template<typename NumericT, int Options>
219 inline vcl_size_t size2(Eigen::Map< Eigen::Matrix<NumericT, Eigen::Dynamic, Eigen::Dynamic, Options> > const & m) { return m.cols(); }
220 template<typename T, int options>
221 inline vcl_size_t size2(Eigen::SparseMatrix<T, options> & m) { return m.cols(); }
222 #endif
223 
224 #ifdef VIENNACL_WITH_MTL4
225 template<typename NumericT, typename T>
226 vcl_size_t size2(mtl::dense2D<NumericT, T> const & m) { return static_cast<vcl_size_t>(m.num_cols()); }
227 template<typename NumericT>
228 vcl_size_t size2(mtl::compressed2D<NumericT> const & m) { return static_cast<vcl_size_t>(m.num_cols()); }
229 #endif
230 
234 //
235 // size: Returns the length of vectors
236 //
238 template<typename VectorType>
239 vcl_size_t size(VectorType const & vec)
240 {
241  return vec.size();
242 }
243 
245 template<typename SparseMatrixType, typename VectorType>
247 {
248  return size1(proxy.lhs());
249 }
250 
251 template<typename T, unsigned int A, typename VectorType>
252 vcl_size_t size(vector_expression<const circulant_matrix<T, A>, const VectorType, op_prod> const & proxy) { return proxy.lhs().size1(); }
253 
254 template<typename T, unsigned int A, typename VectorType>
255 vcl_size_t size(vector_expression<const hankel_matrix<T, A>, const VectorType, op_prod> const & proxy) { return proxy.lhs().size1(); }
256 
257 template<typename T, unsigned int A, typename VectorType>
258 vcl_size_t size(vector_expression<const toeplitz_matrix<T, A>, const VectorType, op_prod> const & proxy) { return proxy.lhs().size1(); }
259 
260 template<typename T, unsigned int A, typename VectorType>
261 vcl_size_t size(vector_expression<const vandermonde_matrix<T, A>, const VectorType, op_prod> const & proxy) { return proxy.lhs().size1(); }
262 
263 template<typename NumericT>
264 vcl_size_t size(vector_expression<const matrix_base<NumericT>, const vector_base<NumericT>, op_prod> const & proxy) //matrix-vector product
265 {
266  return proxy.lhs().size1();
267 }
268 
269 template<typename NumericT, typename LhsT, typename RhsT, typename OpT>
270 vcl_size_t size(vector_expression<const matrix_base<NumericT>, const vector_expression<LhsT, RhsT, OpT>, op_prod> const & proxy) //matrix-vector product
271 {
272  return proxy.lhs().size1();
273 }
274 
275 template<typename NumericT>
276 vcl_size_t size(vector_expression<const matrix_expression<const matrix_base<NumericT>, const matrix_base<NumericT>, op_trans>,
277  const vector_base<NumericT>,
278  op_prod> const & proxy) //transposed matrix-vector product
279 {
280  return proxy.lhs().lhs().size2();
281 }
282 
283 
284 #ifdef VIENNACL_WITH_MTL4
285 template<typename ScalarType>
286 vcl_size_t size(mtl::dense_vector<ScalarType> const & vec) { return vec.used_memory(); }
287 #endif
288 
289 #ifdef VIENNACL_WITH_ARMADILLO
290 template<typename NumericT>
291 inline vcl_size_t size(arma::Mat<NumericT> const & A) { return A.n_elem; }
292 #endif
293 
294 #ifdef VIENNACL_WITH_EIGEN
295 inline vcl_size_t size(Eigen::VectorXf const & v) { return v.rows(); }
296 inline vcl_size_t size(Eigen::VectorXd const & v) { return v.rows(); }
297 #endif
298 
299 template<typename LHS, typename RHS, typename OP>
300 vcl_size_t size(vector_expression<LHS, RHS, OP> const & proxy)
301 {
302  return size(proxy.lhs());
303 }
304 
305 template<typename LHS, typename RHS>
306 vcl_size_t size(vector_expression<LHS, const vector_tuple<RHS>, op_inner_prod> const & proxy)
307 {
308  return proxy.rhs().const_size();
309 }
310 
311 template<typename LhsT, typename RhsT, typename OpT, typename VectorT>
312 vcl_size_t size(vector_expression<const matrix_expression<const LhsT, const RhsT, OpT>,
313  VectorT,
314  op_prod> const & proxy)
315 {
316  return size1(proxy.lhs());
317 }
318 
319 template<typename LhsT, typename RhsT, typename OpT, typename NumericT>
320 vcl_size_t size(vector_expression<const matrix_expression<const LhsT, const RhsT, OpT>,
321  const vector_base<NumericT>,
322  op_prod> const & proxy)
323 {
324  return size1(proxy.lhs());
325 }
326 
327 template<typename LhsT1, typename RhsT1, typename OpT1,
328  typename LhsT2, typename RhsT2, typename OpT2>
329 vcl_size_t size(vector_expression<const matrix_expression<const LhsT1, const RhsT1, OpT1>,
330  const vector_expression<const LhsT2, const RhsT2, OpT2>,
331  op_prod> const & proxy)
332 {
333  return size1(proxy.lhs());
334 }
335 
336 template<typename NumericT>
337 vcl_size_t size(vector_expression<const matrix_base<NumericT>,
338  const matrix_base<NumericT>,
339  op_row_sum> const & proxy)
340 {
341  return size1(proxy.lhs());
342 }
343 
344 template<typename LhsT, typename RhsT, typename OpT>
345 vcl_size_t size(vector_expression<const matrix_expression<const LhsT, const RhsT, OpT>,
346  const matrix_expression<const LhsT, const RhsT, OpT>,
347  op_row_sum> const & proxy)
348 {
349  return size1(proxy.lhs());
350 }
351 
352 template<typename NumericT>
353 vcl_size_t size(vector_expression<const matrix_base<NumericT>,
354  const matrix_base<NumericT>,
355  op_col_sum> const & proxy)
356 {
357  return size2(proxy.lhs());
358 }
359 
360 template<typename LhsT, typename RhsT, typename OpT>
361 vcl_size_t size(vector_expression<const matrix_expression<const LhsT, const RhsT, OpT>,
362  const matrix_expression<const LhsT, const RhsT, OpT>,
363  op_col_sum> const & proxy)
364 {
365  return size2(proxy.lhs());
366 }
367 
370 //
371 // internal_size: Returns the internal (padded) length of vectors
372 //
374 template<typename NumericT>
376 {
377  return vec.internal_size();
378 }
379 
380 
381 //
382 // internal_size1: No. of internal (padded) rows for matrices
383 //
385 template<typename NumericT>
387 
388 
389 //
390 // internal_size2: No. of internal (padded) columns for matrices
391 //
393 template<typename NumericT>
395 
397 template<typename NumericT>
399 {
400  if (mat.row_major())
401  return mat.internal_size2();
402  return mat.internal_size1();
403 }
404 
405 template<typename NumericT>
407 {
408  if (mat.row_major())
409  return mat.stride2();
410  return mat.stride1();
411 }
412 
413 template<typename LHS>
415 {
416  int k = proxy.rhs();
417  int A_size1 = static_cast<int>(size1(proxy.lhs()));
418  int A_size2 = static_cast<int>(size2(proxy.lhs()));
419 
420  int row_depth = std::min(A_size1, A_size1 + k);
421  int col_depth = std::min(A_size2, A_size2 - k);
422 
423  return vcl_size_t(std::min(row_depth, col_depth));
424 }
425 
426 template<typename LHS>
428 {
429  return size2(proxy.lhs());
430 }
431 
432 template<typename LHS>
434 {
435  return size1(proxy.lhs());
436 }
437 
438 } //namespace traits
439 } //namespace viennacl
440 
441 
442 #endif
vcl_size_t nld(matrix_base< NumericT > const &mat)
Definition: size.hpp:406
vcl_size_t internal_size1(matrix_base< NumericT > const &mat)
Helper routine for obtaining the internal number of entries per row of a ViennaCL matrix...
Definition: size.hpp:386
vcl_size_t size1(MatrixType const &mat)
Generic routine for obtaining the number of rows of a matrix (ViennaCL, uBLAS, etc.)
Definition: size.hpp:163
vcl_size_t internal_size2(matrix_base< NumericT > const &mat)
Helper routine for obtaining the internal number of entries per column of a ViennaCL matrix...
Definition: size.hpp:394
size_type stride2() const
Returns the number of columns.
Definition: matrix_def.hpp:234
This file provides the forward declarations for the main types used within ViennaCL.
A dense matrix class.
Definition: forwards.h:375
vcl_size_t internal_size(vector_base< NumericT > const &vec)
Helper routine for obtaining the buffer length of a ViennaCL vector.
Definition: size.hpp:375
lhs_reference_type lhs() const
Get left hand side operand.
Definition: vector.hpp:74
An expression template class that represents a binary operation that yields a vector.
Definition: forwards.h:239
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:201
vcl_size_t size(VectorType const &vec)
Generic routine for obtaining the size of a vector (ViennaCL, uBLAS, etc.)
Definition: size.hpp:239
vcl_size_t ld(matrix_base< NumericT > const &mat)
Helper routine for obtaining the internal number of entries per row of a ViennaCL matrix...
Definition: size.hpp:398
void resize(MatrixType &matrix, vcl_size_t rows, vcl_size_t cols)
Generic resize routine for resizing a matrix (ViennaCL, uBLAS, etc.) to a new size/dimension.
Definition: size.hpp:63
rhs_reference_type rhs() const
Get right hand side operand.
Definition: vector.hpp:77
size_type stride1() const
Returns the number of rows.
Definition: matrix_def.hpp:232
std::size_t vcl_size_t
Definition: forwards.h:75
All the predicates used within ViennaCL. Checks for expressions to be vectors, etc.
bool row_major() const
Definition: matrix_def.hpp:248
size_type internal_size2() const
Returns the internal number of columns. Usually required for launching OpenCL kernels only...
Definition: matrix_def.hpp:240
size_type internal_size1() const
Returns the internal number of rows. Usually required for launching OpenCL kernels only...
Definition: matrix_def.hpp:238
T min(const T &lhs, const T &rhs)
Minimum.
Definition: util.hpp:45
size_type internal_size() const
Returns the internal length of the vector, which is given by size() plus the extra memory due to padd...
Definition: vector_def.hpp:120
A collection of compile time type deductions.