1 #ifndef VIENNACL_TRAITS_SIZE_HPP_
2 #define VIENNACL_TRAITS_SIZE_HPP_
32 #ifdef VIENNACL_WITH_UBLAS
33 #include <boost/numeric/ublas/matrix_sparse.hpp>
34 #include <boost/numeric/ublas/matrix.hpp>
37 #ifdef VIENNACL_WITH_ARMADILLO
41 #ifdef VIENNACL_WITH_EIGEN
43 #include <Eigen/Sparse>
46 #ifdef VIENNACL_WITH_MTL4
47 #include <boost/numeric/mtl/mtl.hpp>
62 template<
typename MatrixType>
65 matrix.resize(rows, cols);
69 template<
typename VectorType>
76 #ifdef VIENNACL_WITH_UBLAS
78 template<
typename ScalarType>
79 void resize(boost::numeric::ublas::compressed_matrix<ScalarType> &
matrix,
83 matrix.resize(rows, cols,
false);
88 #ifdef VIENNACL_WITH_MTL4
89 template<
typename ScalarType>
90 void resize(mtl::compressed2D<ScalarType> & matrix,
94 matrix.change_dim(rows, cols);
97 template<
typename ScalarType>
98 void resize(mtl::dense_vector<ScalarType> & vec,
101 vec.change_dim(new_size);
105 #ifdef VIENNACL_WITH_ARMADILLO
106 template<
typename NumericT>
107 inline void resize(arma::Mat<NumericT> & A,
111 A.resize(new_rows, new_cols);
114 template<
typename NumericT>
115 inline void resize(arma::SpMat<NumericT> & A,
119 A.set_size(new_rows, new_cols);
123 #ifdef VIENNACL_WITH_EIGEN
124 template<
typename NumericT,
int Options>
125 inline void resize(Eigen::Matrix<NumericT, Eigen::Dynamic, Eigen::Dynamic, Options> & m,
129 m.resize(new_rows, new_cols);
132 template<
typename T,
int options>
133 inline void resize(Eigen::SparseMatrix<T, options> & m,
137 m.resize(new_rows, new_cols);
140 inline void resize(Eigen::VectorXf & v,
146 inline void resize(Eigen::VectorXd & v,
161 template<
typename MatrixType>
163 size1(MatrixType
const & mat) {
return mat.size1(); }
166 template<
typename RowType>
168 size1(std::vector< RowType >
const & mat) {
return mat.size(); }
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; }
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>
186 #ifdef VIENNACL_WITH_MTL4
187 template<
typename NumericT,
typename T>
189 template<
typename NumericT>
199 template<
typename MatrixType>
200 typename result_of::size_type<MatrixType>::type
201 size2(MatrixType
const & mat) {
return mat.size2(); }
204 #ifdef VIENNACL_WITH_ARMADILLO
205 template<
typename NumericT>
206 inline vcl_size_t size2(arma::Mat<NumericT>
const & A) {
return A.n_cols; }
207 template<
typename NumericT>
208 inline vcl_size_t size2(arma::SpMat<NumericT>
const & A) {
return A.n_cols; }
211 #ifdef VIENNACL_WITH_EIGEN
212 template<
typename NumericT,
int Options>
213 inline vcl_size_t size2(Eigen::Matrix<NumericT, Eigen::Dynamic, Eigen::Dynamic, Options>
const & m) {
return m.cols(); }
214 template<
typename NumericT,
int Options>
215 inline vcl_size_t size2(Eigen::Map< Eigen::Matrix<NumericT, Eigen::Dynamic, Eigen::Dynamic, Options> >
const & m) {
return m.cols(); }
216 template<
typename T,
int options>
217 inline vcl_size_t size2(Eigen::SparseMatrix<T, options> & m) {
return m.cols(); }
220 #ifdef VIENNACL_WITH_MTL4
221 template<
typename NumericT,
typename T>
223 template<
typename NumericT>
234 template<
typename VectorType>
241 template<
typename SparseMatrixType,
typename VectorType>
247 template<
typename T,
unsigned int A,
typename VectorType>
248 vcl_size_t size(vector_expression<
const circulant_matrix<T, A>,
const VectorType, op_prod>
const & proxy) {
return proxy.
lhs().size1(); }
250 template<
typename T,
unsigned int A,
typename VectorType>
251 vcl_size_t size(vector_expression<
const hankel_matrix<T, A>,
const VectorType, op_prod>
const & proxy) {
return proxy.
lhs().size1(); }
253 template<
typename T,
unsigned int A,
typename VectorType>
254 vcl_size_t size(vector_expression<
const toeplitz_matrix<T, A>,
const VectorType, op_prod>
const & proxy) {
return proxy.
lhs().size1(); }
256 template<
typename T,
unsigned int A,
typename VectorType>
257 vcl_size_t size(vector_expression<
const vandermonde_matrix<T, A>,
const VectorType, op_prod>
const & proxy) {
return proxy.
lhs().size1(); }
259 template<
typename NumericT>
260 vcl_size_t size(vector_expression<
const matrix_base<NumericT>,
const vector_base<NumericT>, op_prod>
const & proxy)
262 return proxy.
lhs().size1();
265 template<
typename NumericT,
typename LhsT,
typename RhsT,
typename OpT>
266 vcl_size_t size(vector_expression<
const matrix_base<NumericT>,
const vector_expression<LhsT, RhsT, OpT>, op_prod>
const & proxy)
268 return proxy.
lhs().size1();
271 template<
typename NumericT>
272 vcl_size_t size(vector_expression<
const matrix_expression<
const matrix_base<NumericT>,
const matrix_base<NumericT>, op_trans>,
273 const vector_base<NumericT>,
274 op_prod>
const & proxy)
276 return proxy.
lhs().lhs().size2();
280 #ifdef VIENNACL_WITH_MTL4
281 template<
typename ScalarType>
282 vcl_size_t size(mtl::dense_vector<ScalarType>
const & vec) {
return vec.used_memory(); }
285 #ifdef VIENNACL_WITH_ARMADILLO
286 template<
typename NumericT>
287 inline vcl_size_t size(arma::Mat<NumericT>
const & A) {
return A.n_elem; }
290 #ifdef VIENNACL_WITH_EIGEN
291 inline vcl_size_t size(Eigen::VectorXf
const & v) {
return v.rows(); }
292 inline vcl_size_t size(Eigen::VectorXd
const & v) {
return v.rows(); }
295 template<
typename LHS,
typename RHS,
typename OP>
298 return size(proxy.lhs());
301 template<
typename LHS,
typename RHS>
302 vcl_size_t size(vector_expression<LHS,
const vector_tuple<RHS>, op_inner_prod>
const & proxy)
304 return proxy.rhs().const_size();
307 template<
typename LhsT,
typename RhsT,
typename OpT,
typename VectorT>
308 vcl_size_t size(vector_expression<
const matrix_expression<const LhsT, const RhsT, OpT>,
310 op_prod>
const & proxy)
312 return size1(proxy.lhs());
315 template<
typename LhsT,
typename RhsT,
typename OpT,
typename NumericT>
316 vcl_size_t size(vector_expression<
const matrix_expression<const LhsT, const RhsT, OpT>,
317 const vector_base<NumericT>,
318 op_prod>
const & proxy)
320 return size1(proxy.lhs());
323 template<
typename LhsT1,
typename RhsT1,
typename OpT1,
324 typename LhsT2,
typename RhsT2,
typename OpT2>
325 vcl_size_t size(vector_expression<
const matrix_expression<const LhsT1, const RhsT1, OpT1>,
326 const vector_expression<const LhsT2, const RhsT2, OpT2>,
327 op_prod>
const & proxy)
329 return size1(proxy.lhs());
332 template<
typename NumericT>
334 const matrix_base<NumericT>,
335 op_row_sum>
const & proxy)
337 return size1(proxy.lhs());
340 template<
typename LhsT,
typename RhsT,
typename OpT>
341 vcl_size_t size(vector_expression<
const matrix_expression<const LhsT, const RhsT, OpT>,
342 const matrix_expression<const LhsT, const RhsT, OpT>,
343 op_row_sum>
const & proxy)
345 return size1(proxy.lhs());
348 template<
typename NumericT>
350 const matrix_base<NumericT>,
351 op_col_sum>
const & proxy)
353 return size2(proxy.lhs());
356 template<
typename LhsT,
typename RhsT,
typename OpT>
357 vcl_size_t size(vector_expression<
const matrix_expression<const LhsT, const RhsT, OpT>,
358 const matrix_expression<const LhsT, const RhsT, OpT>,
359 op_col_sum>
const & proxy)
361 return size2(proxy.lhs());
370 template<
typename NumericT>
381 template<
typename NumericT>
389 template<
typename NumericT>
393 template<
typename NumericT>
401 template<
typename NumericT>
409 template<
typename LHS>
413 int A_size1 =
static_cast<int>(
size1(proxy.
lhs()));
414 int A_size2 =
static_cast<int>(
size2(proxy.
lhs()));
416 int row_depth =
std::min(A_size1, A_size1 + k);
417 int col_depth =
std::min(A_size2, A_size2 - k);
422 template<
typename LHS>
428 template<
typename LHS>
vcl_size_t nld(matrix_base< NumericT > const &mat)
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...
vcl_size_t size1(MatrixType const &mat)
Generic routine for obtaining the number of rows of a matrix (ViennaCL, uBLAS, etc.)
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...
size_type stride2() const
Returns the number of columns.
This file provides the forward declarations for the main types used within ViennaCL.
vcl_size_t internal_size(vector_base< NumericT > const &vec)
Helper routine for obtaining the buffer length of a ViennaCL vector.
lhs_reference_type lhs() const
Get left hand side operand.
An expression template class that represents a binary operation that yields a vector.
result_of::size_type< MatrixType >::type size2(MatrixType const &mat)
Generic routine for obtaining the number of columns of a matrix (ViennaCL, uBLAS, etc...
Main namespace in ViennaCL. Holds all the basic types such as vector, matrix, etc. and defines operations upon them.
vcl_size_t size(VectorType const &vec)
Generic routine for obtaining the size of a vector (ViennaCL, uBLAS, etc.)
vcl_size_t ld(matrix_base< NumericT > const &mat)
Helper routine for obtaining the internal number of entries per row of a ViennaCL matrix...
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.
rhs_reference_type rhs() const
Get right hand side operand.
size_type stride1() const
Returns the number of rows.
All the predicates used within ViennaCL. Checks for expressions to be vectors, etc.
size_type internal_size2() const
Returns the internal number of columns. Usually required for launching OpenCL kernels only...
size_type internal_size1() const
Returns the internal number of rows. Usually required for launching OpenCL kernels only...
T min(const T &lhs, const T &rhs)
Minimum.
size_type internal_size() const
Returns the internal length of the vector, which is given by size() plus the extra memory due to padd...
A collection of compile time type deductions.