1 #ifndef VIENNACL_MATRIX_HPP_
2 #define VIENNACL_MATRIX_HPP_
53 template<
typename LHS,
typename RHS,
typename OP>
54 class matrix_expression
56 typedef typename viennacl::result_of::reference_if_nonscalar<LHS>::type lhs_reference_type;
57 typedef typename viennacl::result_of::reference_if_nonscalar<RHS>::type rhs_reference_type;
66 LHS &
lhs()
const {
return lhs_; }
69 RHS &
rhs()
const {
return rhs_; }
77 lhs_reference_type lhs_;
79 rhs_reference_type rhs_;
91 template<
typename ROWCOL,
typename MatrixT>
100 vcl_size_t start_col) : mat_(mat), row_(start_row), col_(start_col) {}
104 self_type
operator++(
int) { self_type tmp = *
this; ++(*this);
return tmp; }
106 bool operator==(self_type
const & other) {
return (row_ == other.row_) && (col_ == other.col_); }
107 bool operator!=(self_type
const & other) {
return !(*
this == other); }
126 template<
class NumericT,
typename SizeT,
typename DistanceT>
128 : size1_(rows), size2_(columns), start1_(0), start2_(0), stride1_(1), stride2_(1),
131 row_major_fixed_(true), row_major_(is_row_major)
133 if (rows > 0 && columns > 0)
142 template<
class NumericT,
typename SizeT,
typename DistanceT>
143 template<
typename LHS,
typename RHS,
typename OP>
145 size1_(
viennacl::traits::
size1(proxy)), size2_(
viennacl::traits::
size2(proxy)), start1_(0), start2_(0), stride1_(1), stride2_(1),
160 template<
class NumericT,
typename SizeT,
typename DistanceT>
165 : size1_(mat_size1), size2_(mat_size2),
166 start1_(mat_start1), start2_(mat_start2),
167 stride1_(mat_stride1), stride2_(mat_stride2),
168 internal_size1_(mat_internal_size1), internal_size2_(mat_internal_size2),
169 row_major_fixed_(true), row_major_(is_row_major)
173 #ifdef VIENNACL_WITH_CUDA
175 elements_.cuda_handle().reset(reinterpret_cast<char*>(ptr_to_mem));
176 elements_.cuda_handle().inc();
191 #ifdef VIENNACL_WITH_OPENCL
192 template<
class NumericT,
typename SizeT,
typename DistanceT>
194 : size1_(rows), size2_(columns),
195 start1_(0), start2_(0),
196 stride1_(1), stride2_(1),
197 internal_size1_(rows), internal_size2_(columns),
198 row_major_fixed_(true), row_major_(is_row_major)
201 elements_.opencl_handle() = mem;
202 elements_.opencl_handle().inc();
203 elements_.opencl_handle().context(ctx.opencl_context());
207 template<
class NumericT,
typename SizeT,
typename DistanceT>
209 size_type mat_size1, size_type mat_start1, size_type mat_stride1, size_type mat_internal_size1,
210 size_type mat_size2, size_type mat_start2, size_type mat_stride2, size_type mat_internal_size2,
212 : size1_(mat_size1), size2_(mat_size2),
213 start1_(mat_start1), start2_(mat_start2),
214 stride1_(mat_stride1), stride2_(mat_stride2),
215 internal_size1_(mat_internal_size1), internal_size2_(mat_internal_size2),
216 row_major_fixed_(true), row_major_(is_row_major)
219 elements_.opencl_handle() = mem;
220 elements_.opencl_handle().inc();
221 elements_.opencl_handle().context(ctx.opencl_context());
227 template<
class NumericT,
typename SizeT,
typename DistanceT>
229 size1_(other.
size1()), size2_(other.
size2()), start1_(0), start2_(0), stride1_(1), stride2_(1),
232 row_major_fixed_(true), row_major_(other.
row_major())
244 template<
typename NumericT,
typename SizeT,
typename DistanceT>
245 template<
typename OtherNumericT>
247 size1_(other.
size1()), size2_(other.
size2()), start1_(0), start2_(0), stride1_(1), stride2_(1),
250 row_major_fixed_(true), row_major_(other.
row_major())
261 template<
class NumericT,
typename SizeT,
typename DistanceT>
271 if (!row_major_fixed_)
282 template<
class NumericT,
typename SizeT,
typename DistanceT>
283 template<
typename OtherNumericT>
290 if (!row_major_fixed_)
303 template<
class NumericT,
typename SizeT,
typename DistanceT>
304 template<
typename LHS,
typename RHS,
typename OP>
309 &&
bool(
"Incompatible matrix sizes!"));
314 internal_size1_ = viennacl::tools::align_to_multiple<size_type>(size1_, dense_padding_size);
315 internal_size2_ = viennacl::tools::align_to_multiple<size_type>(size2_, dense_padding_size);
316 if (!row_major_fixed_)
319 if (size1_ != internal_size1_ || size2_ != internal_size2_)
331 template<
class NumericT,
typename SizeT,
typename DistanceT>
338 internal_size1_ = viennacl::tools::align_to_multiple<size_type>(size1_, dense_padding_size);
339 internal_size2_ = viennacl::tools::align_to_multiple<size_type>(size2_, dense_padding_size);
340 if (!row_major_fixed_)
344 if (
handle() == proxy.lhs().handle() )
348 if ( proxy.lhs().size1() != proxy.lhs().size2() )
349 this->
resize(proxy.lhs().size2(), proxy.lhs().size1());
350 elements_ = temp.handle();
354 if ( proxy.lhs().size1() != proxy.lhs().size2() )
355 this->
resize(proxy.lhs().size2(), proxy.lhs().size1());
361 template<
class NumericT,
typename SizeT,
typename DistanceT>
362 template<
typename LHS,
typename RHS,
typename OP>
367 &&
bool(
"Incompatible matrix sizes!"));
368 assert( (
size1() > 0) &&
bool(
"Vector not yet initialized!") );
369 assert( (
size2() > 0) &&
bool(
"Vector not yet initialized!") );
376 template<
class NumericT,
typename SizeT,
typename DistanceT>
377 template<
typename LHS,
typename RHS,
typename OP>
382 &&
bool(
"Incompatible matrix sizes!"));
383 assert( (
size1() > 0) &&
bool(
"Vector not yet initialized!") );
384 assert( (
size2() > 0) &&
bool(
"Vector not yet initialized!") );
392 template<
class NumericT,
typename SizeT,
typename DistanceT>
395 assert( (m.
size1() == size1_ || size1_ == 0) &&
bool(
"Size mismatch!") );
396 assert( (m.
size2() == size2_ || size2_ == 0) &&
bool(
"Size mismatch!") );
402 internal_size1_ = viennacl::tools::align_to_multiple<size_type>(size1_, dense_padding_size);
403 internal_size2_ = viennacl::tools::align_to_multiple<size_type>(size2_, dense_padding_size);
420 template<
class NumericT,
typename SizeT,
typename DistanceT>
423 assert( (m.
size1() == size1_ || size1_ == 0) &&
bool(
"Size mismatch!") );
424 assert( (m.
size2() == size2_ || size2_ == 0) &&
bool(
"Size mismatch!") );
430 internal_size1_ = viennacl::tools::align_to_multiple<size_type>(size1_, dense_padding_size);
431 internal_size2_ = viennacl::tools::align_to_multiple<size_type>(size2_, dense_padding_size);
445 template<
class NumericT,
typename SizeT,
typename DistanceT>
448 assert( (m.
size1() == size1_ || size1_ == 0) &&
bool(
"Size mismatch!") );
449 assert( (m.
size2() == size2_ || size2_ == 0) &&
bool(
"Size mismatch!") );
455 internal_size1_ = viennacl::tools::align_to_multiple<size_type>(size1_, dense_padding_size);
456 internal_size2_ = viennacl::tools::align_to_multiple<size_type>(size2_, dense_padding_size);
476 template<
class NumericT,
typename SizeT,
typename DistanceT>
486 template<
class NumericT,
typename SizeT,
typename DistanceT>
497 template<
class NumericT,
typename SizeT,
typename DistanceT>
501 *
this,
NumericT(1.0), 1,
false,
false,
502 other,
NumericT(1.0), 1,
false,
false);
506 template<
class NumericT,
typename SizeT,
typename DistanceT>
510 *
this,
NumericT(1.0), 1,
false,
false,
511 other,
NumericT(1.0), 1,
false,
true);
516 template<
class NumericT,
typename SizeT,
typename DistanceT>
520 *
this,
NumericT(val), 1,
false,
false);
525 template<
class NumericT,
typename SizeT,
typename DistanceT>
529 *
this,
NumericT(val), 1,
false,
false);
534 template<
class NumericT,
typename SizeT,
typename DistanceT>
538 *
this,
NumericT(val), 1,
false,
false);
543 template<
class NumericT,
typename SizeT,
typename DistanceT>
547 *
this,
NumericT(val), 1,
false,
false);
552 template<
class NumericT,
typename SizeT,
typename DistanceT>
556 *
this,
NumericT(val), 1,
false,
false);
561 template<
class NumericT,
typename SizeT,
typename DistanceT>
565 *
this,
NumericT(val), 1,
false,
false);
572 template<
class NumericT,
typename SizeT,
typename DistanceT>
576 *
this,
NumericT(val), 1,
true,
false);
581 template<
class NumericT,
typename SizeT,
typename DistanceT>
585 *
this,
NumericT(val), 1,
true,
false);
590 template<
class NumericT,
typename SizeT,
typename DistanceT>
594 *
this,
NumericT(val), 1,
true,
false);
599 template<
class NumericT,
typename SizeT,
typename DistanceT>
603 *
this,
NumericT(val), 1,
true,
false);
608 template<
class NumericT,
typename SizeT,
typename DistanceT>
612 *
this,
NumericT(val), 1,
true,
false);
617 template<
class NumericT,
typename SizeT,
typename DistanceT>
621 *
this,
NumericT(val), 1,
true,
false);
627 template<
class NumericT,
typename SizeT,
typename DistanceT>
633 template<
class NumericT,
typename SizeT,
typename DistanceT>
637 template<
class NumericT,
typename SizeT,
typename DistanceT>
640 assert( (rows > 0 && columns > 0) &&
bool(
"Check failed in matrix::resize(): Number of rows and columns must be positive!"));
649 std::vector< NumericT > new_entries( viennacl::tools::align_to_multiple<vcl_size_t>(rows, dense_padding_size)
650 * viennacl::tools::align_to_multiple<vcl_size_t>(columns, dense_padding_size));
661 new_entries[
row_major::mem_index(i, j, viennacl::tools::align_to_multiple<vcl_size_t>(rows, dense_padding_size), viennacl::tools::align_to_multiple<vcl_size_t>(columns, dense_padding_size))]
664 new_entries[
column_major::mem_index(i, j, viennacl::tools::align_to_multiple<vcl_size_t>(rows, dense_padding_size), viennacl::tools::align_to_multiple<vcl_size_t>(columns, dense_padding_size))]
672 internal_size1_ = viennacl::tools::align_to_multiple<size_type>(size1_, dense_padding_size);
673 internal_size2_ = viennacl::tools::align_to_multiple<size_type>(size2_, dense_padding_size);
680 internal_size1_ = viennacl::tools::align_to_multiple<size_type>(size1_, dense_padding_size);
681 internal_size2_ = viennacl::tools::align_to_multiple<size_type>(size2_, dense_padding_size);
695 template<
class NumericT,
typename F,
unsigned int AlignmentV>
704 explicit matrix() : base_type(static_cast<bool>(
viennacl::is_row_major<F>::value)) {}
722 : base_type(ptr_to_mem, mem_type,
725 viennacl::is_row_major<F>::value) {}
738 size_type rows, size_type internal_row_count,
739 size_type cols, size_type internal_col_count)
740 : base_type(ptr_to_mem, mem_type,
741 rows, 0, 1, internal_row_count,
742 cols, 0, 1, internal_col_count,
743 true,
viennacl::is_row_major<F>::value) {}
745 #ifdef VIENNACL_WITH_OPENCL
746 explicit matrix(cl_mem mem, size_type rows, size_type columns) : base_type(mem, rows, columns,
viennacl::is_row_major<F>::value) {}
749 template<
typename LHS,
typename RHS,
typename OP>
794 using base_type::operator=;
797 template<
typename OtherNumericT,
typename F2>
800 template<
typename OtherNumericT,
typename F2>
803 template<
typename OtherNumericT,
typename F2>
813 void resize(size_type rows, size_type columns,
bool preserve =
true)
827 template<
class NumericT>
828 std::ostream & operator<<(std::ostream & s, const matrix_base<NumericT> & gpu_matrix)
832 std::vector<NumericT> tmp(gpu_matrix.internal_size());
835 s <<
"[" << gpu_matrix.size1() <<
"," << gpu_matrix.size2() <<
"]";
838 for (size_type i = 0; i < gpu_matrix.size1(); ++i)
841 for (size_type j = 0; j < gpu_matrix.size2(); ++j)
843 if (gpu_matrix.row_major())
844 s << tmp[
row_major::mem_index(i * gpu_matrix.stride1() + gpu_matrix.start1(), j * gpu_matrix.stride2() + gpu_matrix.start2(), gpu_matrix.internal_size1(), gpu_matrix.internal_size2())];
846 s << tmp[
column_major::mem_index(i * gpu_matrix.stride1() + gpu_matrix.start1(), j * gpu_matrix.stride2() + gpu_matrix.start2(), gpu_matrix.internal_size1(), gpu_matrix.internal_size2())];
848 if (j < gpu_matrix.size2() - 1)
852 if (i < gpu_matrix.size1() - 1)
864 template<
typename LHS,
typename RHS,
typename OP>
865 std::ostream & operator<<(std::ostream & s, const matrix_expression<LHS, RHS, OP> & expr)
875 template<
typename NumericT>
876 matrix_expression< const matrix_base<NumericT>,
const matrix_base<NumericT>, op_trans>
877 trans(
const matrix_base<NumericT> & mat)
883 template<
typename NumericT>
884 vector_expression< const matrix_base<NumericT>,
const int, op_matrix_diag>
885 diag(
const matrix_base<NumericT> & A,
int k = 0)
890 template<
typename NumericT>
891 matrix_expression< const vector_base<NumericT>,
const int, op_vector_diag>
898 template<
typename NumericT,
typename F>
899 vector_expression< const matrix_base<NumericT, F>,
const unsigned int, op_row>
906 template<
typename NumericT,
typename F>
907 vector_expression< const matrix_base<NumericT, F>,
const unsigned int, op_column>
923 template<
typename CPUMatrixT,
typename NumericT,
typename F,
unsigned int AlignmentV>
924 void copy(
const CPUMatrixT & cpu_matrix,
931 if (gpu_matrix.
size1() == 0 || gpu_matrix.
size2() == 0)
933 gpu_matrix.
resize(cpu_matrix.size1(),
934 cpu_matrix.size2(),
false);
937 assert( (gpu_matrix.
size1() == cpu_matrix.size1()) && (gpu_matrix.
size2() == cpu_matrix.size2()) &&
bool(
"Matrix dimensions mismatch.") );
940 for (size_type i = 0; i < gpu_matrix.
size1(); ++i)
942 for (size_type j = 0; j < gpu_matrix.
size2(); ++j)
959 template<
typename NumericT,
typename A1,
typename A2,
typename F,
unsigned int AlignmentV>
960 void copy(
const std::vector< std::vector<NumericT, A1>, A2> & cpu_matrix,
965 if (gpu_matrix.
size1() == 0 || gpu_matrix.
size2() == 0)
967 gpu_matrix.
resize(cpu_matrix.size(),
968 cpu_matrix[0].size(),
972 assert( (gpu_matrix.
size1() == cpu_matrix.size()) &&
bool(
"Matrix dimensions mismatch.") );
975 for (size_type i = 0; i < gpu_matrix.
size1(); ++i)
977 assert( (gpu_matrix.
size2() == cpu_matrix[i].size()) &&
bool(
"Matrix dimensions mismatch.") );
979 for (size_type j = 0; j < gpu_matrix.
size2(); ++j)
999 template<
typename NumericT,
typename F,
unsigned int AlignmentV>
1001 NumericT * cpu_matrix_end,
1008 assert( (gpu_matrix.
internal_size() >=
static_cast<vcl_size_t>(cpu_matrix_end - cpu_matrix_begin)) && bool(
"fast_copy(): Matrix not large enough to fit data!"));
1013 #ifdef VIENNACL_WITH_ARMADILLO
1019 template<
typename NumericT,
typename F,
unsigned int AlignmentV>
1020 void copy(arma::Mat<NumericT>
const & arma_matrix,
1025 if (vcl_matrix.
size1() == 0 || vcl_matrix.
size2() == 0)
1027 vcl_matrix.
resize(arma_matrix.n_rows,
1033 assert( (vcl_matrix.
size1() ==
static_cast<vcl_size_t>(arma_matrix.n_rows))
1034 && (vcl_matrix.
size2() ==
static_cast<vcl_size_t>(arma_matrix.n_cols))
1035 && bool(
"matrix size mismatch")
1041 for (size_type j = 0; j < vcl_matrix.
size2(); ++j)
1042 for (size_type i = 0; i < vcl_matrix.
size1(); ++i)
1050 #ifdef VIENNACL_WITH_EIGEN
1053 template<
typename EigenMatrixTypeT,
typename NumericT,
typename F,
unsigned int AlignmentV>
1054 void copy_from_eigen_matrix(EigenMatrixTypeT
const & cpu_matrix,
1059 if (gpu_matrix.
size1() == 0 || gpu_matrix.
size2() == 0)
1061 gpu_matrix.
resize(cpu_matrix.rows(),
1067 assert( (gpu_matrix.
size1() ==
static_cast<vcl_size_t>(cpu_matrix.rows()))
1069 && bool(
"matrix size mismatch")
1074 for (size_type i = 0; i < gpu_matrix.
size1(); ++i)
1076 for (size_type j = 0; j < gpu_matrix.
size2(); ++j)
1090 template<
typename NumericT,
int EigenOptions,
typename F,
unsigned int AlignmentV>
1091 void copy(Eigen::Matrix<NumericT, Eigen::Dynamic, Eigen::Dynamic, EigenOptions>
const & cpu_matrix,
1094 detail::copy_from_eigen_matrix(cpu_matrix, vcl_matrix);
1102 template<
typename NumericT,
int EigenOptions,
int EigenMatTypeV,
typename EigenStr
ideT,
typename F,
unsigned int AlignmentV>
1103 void copy(Eigen::Map<Eigen::Matrix<NumericT, Eigen::Dynamic, Eigen::Dynamic, EigenOptions>, EigenMatTypeV, EigenStrideT>
const & cpu_matrix,
1106 detail::copy_from_eigen_matrix(cpu_matrix, vcl_matrix);
1110 #ifdef VIENNACL_WITH_MTL4
1116 template<
typename NumericT,
typename T,
typename F,
unsigned int AlignmentV>
1117 void copy(
const mtl::dense2D<NumericT, T>& cpu_matrix,
1118 matrix<NumericT, F, AlignmentV> & gpu_matrix)
1120 typedef typename matrix<NumericT, F, AlignmentV>::size_type size_type;
1122 if (gpu_matrix.size1() == 0 || gpu_matrix.size2() == 0)
1124 gpu_matrix.resize(cpu_matrix.num_rows(),
1125 cpu_matrix.num_cols(),
1130 assert( (gpu_matrix.size1() == cpu_matrix.num_rows())
1131 && (gpu_matrix.size2() == cpu_matrix.num_cols())
1132 &&
bool(
"matrix size mismatch")
1136 std::vector<NumericT> data(gpu_matrix.internal_size());
1137 for (size_type i = 0; i < gpu_matrix.size1(); ++i)
1139 for (size_type j = 0; j < gpu_matrix.size2(); ++j)
1140 data[F::mem_index(i, j, gpu_matrix.internal_size1(), gpu_matrix.internal_size2())] = cpu_matrix[i][j];
1159 template<
typename CPUMatrixT,
typename NumericT,
typename F,
unsigned int AlignmentV>
1161 CPUMatrixT & cpu_matrix )
1165 if ( (gpu_matrix.
size1() > 0) && (gpu_matrix.
size2() > 0) )
1169 std::vector<NumericT> temp_buffer(gpu_matrix.
internal_size());
1173 for (size_type i = 0; i < gpu_matrix.
size1(); ++i)
1176 for (size_type j = 0; j < gpu_matrix.
size2(); ++j)
1188 template<
typename NumericT,
typename A1,
typename A2,
typename F,
unsigned int AlignmentV>
1190 std::vector< std::vector<NumericT, A1>, A2> & cpu_matrix)
1194 if ( (gpu_matrix.
size1() > 0) && (gpu_matrix.
size2() > 0) )
1196 assert( (cpu_matrix.size() == gpu_matrix.
size1()) &&
bool(
"Matrix dimensions mismatch: rows"));
1198 std::vector<NumericT> temp_buffer(gpu_matrix.
internal_size());
1202 for (size_type i = 0; i < gpu_matrix.
size1(); ++i)
1204 assert( (cpu_matrix[i].
size() == gpu_matrix.
size2()) &&
bool(
"Matrix dimensions mismatch: columns"));
1206 for (size_type j = 0; j < gpu_matrix.
size2(); ++j)
1220 template<
typename NumericT,
typename F,
unsigned int AlignmentV>
1222 NumericT * cpu_matrix_begin)
1234 template<
typename LHS1,
typename RHS1,
typename OP1,
1235 typename LHS2,
typename RHS2,
typename OP2>
1236 matrix_expression< const matrix_expression<const LHS1, const RHS1, OP1>,
1237 const matrix_expression<const LHS2, const RHS2, OP2>,
1240 matrix_expression<const LHS2, const RHS2, OP2>
const & proxy2)
1244 &&
bool(
"Incompatible matrix sizes!"));
1245 return matrix_expression< const matrix_expression<const LHS1, const RHS1, OP1>,
1246 const matrix_expression<const LHS2, const RHS2, OP2>,
1250 template<
typename LHS1,
typename RHS1,
typename OP1,
1252 matrix_expression< const matrix_expression<const LHS1, const RHS1, OP1>,
1253 const matrix_base<NumericT>,
1256 matrix_base<NumericT>
const & proxy2)
1260 &&
bool(
"Incompatible matrix sizes!"));
1261 return matrix_expression< const matrix_expression<const LHS1, const RHS1, OP1>,
1262 const matrix_base<NumericT>,
1267 typename LHS2,
typename RHS2,
typename OP2>
1268 matrix_expression< const matrix_base<NumericT>,
1269 const matrix_expression<const LHS2, const RHS2, OP2>,
1272 matrix_expression<const LHS2, const RHS2, OP2>
const & proxy2)
1276 &&
bool(
"Incompatible matrix sizes!"));
1277 return matrix_expression< const matrix_base<NumericT>,
1278 const matrix_expression<const LHS2, const RHS2, OP2>,
1283 template<
typename NumericT>
1284 matrix_expression< const matrix_base<NumericT>,
const matrix_base<NumericT>, op_add >
1285 operator + (
const matrix_base<NumericT> &
m1,
const matrix_base<NumericT> & m2)
1287 return matrix_expression< const matrix_base<NumericT>,
1288 const matrix_base<NumericT>,
1294 template<
typename LHS1,
typename RHS1,
typename OP1,
1295 typename LHS2,
typename RHS2,
typename OP2>
1296 matrix_expression< const matrix_expression<const LHS1, const RHS1, OP1>,
1297 const matrix_expression<const LHS2, const RHS2, OP2>,
1300 matrix_expression<const LHS2, const RHS2, OP2>
const & proxy2)
1304 &&
bool(
"Incompatible matrix sizes!"));
1305 return matrix_expression< const matrix_expression<const LHS1, const RHS1, OP1>,
1306 const matrix_expression<const LHS2, const RHS2, OP2>,
1310 template<
typename LHS1,
typename RHS1,
typename OP1,
1312 matrix_expression< const matrix_expression<const LHS1, const RHS1, OP1>,
1313 const matrix_base<NumericT>,
1316 matrix_base<NumericT>
const & proxy2)
1320 &&
bool(
"Incompatible matrix sizes!"));
1321 return matrix_expression< const matrix_expression<const LHS1, const RHS1, OP1>,
1322 const matrix_base<NumericT>,
1327 typename LHS2,
typename RHS2,
typename OP2>
1328 matrix_expression< const matrix_base<NumericT>,
1329 const matrix_expression<const LHS2, const RHS2, OP2>,
1332 matrix_expression<const LHS2, const RHS2, OP2>
const & proxy2)
1336 &&
bool(
"Incompatible matrix sizes!"));
1337 return matrix_expression< const matrix_base<NumericT>,
1338 const matrix_expression<const LHS2, const RHS2, OP2>,
1343 template<
typename NumericT>
1344 matrix_expression< const matrix_base<NumericT>,
const matrix_base<NumericT>, op_sub >
1345 operator - (
const matrix_base<NumericT> &
m1,
const matrix_base<NumericT> & m2)
1347 return matrix_expression< const matrix_base<NumericT>,
1348 const matrix_base<NumericT>,
1360 template<
typename S1,
typename NumericT>
1362 matrix_expression< const matrix_base<NumericT>,
const S1, op_mult>
1366 return matrix_expression< const matrix_base<NumericT>,
const S1,
op_mult>(
m1, value);
1370 template<
typename NumericT>
1371 matrix_expression< const matrix_base<NumericT>,
const NumericT, op_mult>
1378 template<
typename NumericT>
1379 matrix_expression< const matrix_base<NumericT>,
const NumericT, op_mult>
1386 template<
typename NumericT>
1387 matrix_expression< const matrix_base<NumericT>,
const NumericT, op_mult>
1394 template<
typename NumericT>
1395 matrix_expression< const matrix_base<NumericT>,
const NumericT, op_mult>
1402 template<
typename NumericT>
1403 matrix_expression< const matrix_base<NumericT>,
const NumericT, op_mult>
1410 template<
typename NumericT>
1411 matrix_expression< const matrix_base<NumericT>,
const NumericT, op_mult>
1424 template<
typename LHS,
typename RHS,
typename OP,
typename S1>
1426 matrix_expression< const matrix_expression< LHS, RHS, OP>,
const S1, op_mult> >::type
1430 return matrix_expression< const matrix_expression< LHS, RHS, OP>,
const S1,
op_mult>(proxy, val);
1439 template<
typename S1,
typename LHS,
typename RHS,
typename OP>
1441 matrix_expression< const matrix_expression< LHS, RHS, OP>,
const S1, op_mult> >::type
1445 return matrix_expression< const matrix_expression< LHS, RHS, OP>,
const S1,
op_mult>(proxy, val);
1450 template<
typename NumericT,
typename S1>
1452 matrix_expression< const matrix_base<NumericT>,
const S1, op_mult> >::type
1455 return matrix_expression< const matrix_base<NumericT>,
const S1,
op_mult>(
m1,
s1);
1459 template<
typename NumericT>
1460 matrix_expression< const matrix_base<NumericT>,
const NumericT, op_mult>
1467 template<
typename NumericT>
1468 matrix_expression< const matrix_base<NumericT>,
const NumericT, op_mult>
1475 template<
typename NumericT>
1476 matrix_expression< const matrix_base<NumericT>,
const NumericT, op_mult>
1483 template<
typename NumericT>
1484 matrix_expression< const matrix_base<NumericT>,
const NumericT, op_mult>
1491 template<
typename NumericT>
1492 matrix_expression< const matrix_base<NumericT>,
const NumericT, op_mult>
1499 template<
typename NumericT>
1500 matrix_expression< const matrix_base<NumericT>,
const NumericT, op_mult>
1510 template<
typename NumericT,
typename S1>
1516 m1, gpu_val, 1,
false, is_sign_flip ?
true :
false);
1521 template<
typename NumericT>
1522 matrix_base<NumericT> &
1526 m1,
NumericT(gpu_val), 1,
false,
false);
1531 template<
typename NumericT>
1532 matrix_base<NumericT> &
1536 m1,
NumericT(gpu_val), 1,
false,
false);
1541 template<
typename NumericT>
1542 matrix_base<NumericT> &
1546 m1,
NumericT(gpu_val), 1,
false,
false);
1551 template<
typename NumericT>
1552 matrix_base<NumericT> &
1556 m1,
NumericT(gpu_val), 1,
false,
false);
1561 template<
typename NumericT>
1562 matrix_base<NumericT> &
1566 m1,
NumericT(gpu_val), 1,
false,
false);
1571 template<
typename NumericT>
1572 matrix_base<NumericT> &
1576 m1,
NumericT(gpu_val), 1,
false,
false);
1590 template<
typename LHS,
typename RHS,
typename OP,
typename S1>
1592 matrix_expression< const matrix_expression<const LHS, const RHS, OP>,
const S1, op_div> >::type
1596 return matrix_expression< const matrix_expression<const LHS, const RHS, OP>,
const S1,
op_div>(proxy, val);
1601 template<
typename NumericT,
typename S1>
1603 matrix_expression< const matrix_base<NumericT>,
const S1, op_div> >::type
1606 return matrix_expression< const matrix_base<NumericT>,
const S1,
op_div>(
m1,
s1);
1610 template<
typename NumericT>
1611 matrix_expression< const matrix_base<NumericT>,
const NumericT, op_div>
1618 template<
typename NumericT>
1619 matrix_expression< const matrix_base<NumericT>,
const NumericT, op_div>
1626 template<
typename NumericT>
1627 matrix_expression< const matrix_base<NumericT>,
const NumericT, op_div>
1634 template<
typename NumericT>
1635 matrix_expression< const matrix_base<NumericT>,
const NumericT, op_div>
1642 template<
typename NumericT>
1643 matrix_expression< const matrix_base<NumericT>,
const NumericT, op_div>
1650 template<
typename NumericT>
1651 matrix_expression< const matrix_base<NumericT>,
const NumericT, op_div>
1662 template<
typename NumericT,
typename S1>
1667 m1, gpu_val, 1,
true,
false);
1672 template<
typename NumericT>
1673 matrix_base<NumericT> &
1677 m1,
NumericT(gpu_val), 1,
true,
false);
1682 template<
typename NumericT>
1683 matrix_base<NumericT> &
1687 m1, gpu_val, 1,
true,
false);
1692 template<
typename NumericT>
1693 matrix_base<NumericT> &
1697 m1, gpu_val, 1,
true,
false);
1702 template<
typename NumericT>
1703 matrix_base<NumericT> &
1707 m1, gpu_val, 1,
true,
false);
1712 template<
typename NumericT>
1713 matrix_base<NumericT> &
1717 m1, gpu_val, 1,
true,
false);
1722 template<
typename NumericT>
1723 matrix_base<NumericT> &
1727 m1, gpu_val, 1,
true,
false);
1736 template<
typename NumericT,
typename S1>
1745 return viennacl::matrix_expression< const viennacl::matrix_expression< const vector_base<NumericT>,
const vector_base<NumericT>, op_prod>,
1750 template<
typename NumericT,
typename S1>
1752 viennacl::matrix_expression< const viennacl::matrix_expression< const vector_base<NumericT>,
const vector_base<NumericT>, op_prod>,
1759 return viennacl::matrix_expression< const viennacl::matrix_expression< const vector_base<NumericT>,
const vector_base<NumericT>, op_prod>,
1765 template<
typename NumericT,
typename S1>
1767 viennacl::matrix_expression< const viennacl::matrix_expression< const vector_base<NumericT>,
const vector_base<NumericT>, op_prod>,
1774 return viennacl::matrix_expression< const viennacl::matrix_expression< const vector_base<NumericT>,
const vector_base<NumericT>, op_prod>,
1779 template<
typename NumericT,
typename S1>
1781 viennacl::matrix_expression< const viennacl::matrix_expression< const vector_base<NumericT>,
const vector_base<NumericT>, op_prod>,
1788 return viennacl::matrix_expression< const viennacl::matrix_expression< const vector_base<NumericT>,
const vector_base<NumericT>, op_prod>,
1807 template<
typename T>
1808 struct op_executor<matrix_base<T>,
op_assign, matrix_base<T> >
1810 static void apply(matrix_base<T> & lhs, matrix_base<T>
const & rhs)
1817 template<
typename T>
1818 struct op_executor<matrix_base<T>,
op_assign, matrix_expression<const matrix_base<T>, const matrix_base<T>, op_trans> >
1820 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_base<T>,
const matrix_base<T>, op_trans>
const & rhs)
1822 matrix_base<T> temp(rhs);
1829 template<
typename T>
1830 struct op_executor<matrix_base<T>, op_inplace_add, matrix_base<T> >
1832 static void apply(matrix_base<T> & lhs, matrix_base<T>
const & rhs)
1834 viennacl::linalg::ambm(lhs, lhs, T(1), 1,
false,
false, rhs, T(1), 1,
false,
false);
1839 template<
typename T>
1840 struct op_executor<matrix_base<T>, op_inplace_add, matrix_expression<const matrix_base<T>, const matrix_base<T>, op_trans> >
1842 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_base<T>,
const matrix_base<T>, op_trans>
const & rhs)
1844 matrix_base<T> temp(rhs);
1845 viennacl::linalg::ambm(lhs, lhs, T(1), 1,
false,
false, temp, T(1), 1,
false,
false);
1850 template<
typename T>
1851 struct op_executor<matrix_base<T>, op_inplace_sub, matrix_base<T> >
1853 static void apply(matrix_base<T> & lhs, matrix_base<T>
const & rhs)
1855 viennacl::linalg::ambm(lhs, lhs, T(1), 1,
false,
false, rhs, T(1), 1,
false,
true);
1860 template<
typename T>
1861 struct op_executor<matrix_base<T>, op_inplace_sub, matrix_expression<const matrix_base<T>, const matrix_base<T>, op_trans> >
1863 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_base<T>,
const matrix_base<T>, op_trans>
const & rhs)
1865 matrix_base<T> temp(rhs);
1866 viennacl::linalg::ambm(lhs, lhs, T(1), 1,
false,
false, temp, T(1), 1,
false,
true);
1874 template<
typename T,
typename ScalarType>
1875 struct op_executor<matrix_base<T>,
op_assign, matrix_expression<const matrix_base<T>, const
ScalarType, op_mult> >
1877 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_base<T>,
const ScalarType, op_mult>
const & proxy)
1884 template<
typename T,
typename ScalarType>
1885 struct op_executor<matrix_base<T>, op_inplace_add, matrix_expression<const matrix_base<T>, const
ScalarType, op_mult> >
1887 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_base<T>,
const ScalarType, op_mult>
const & proxy)
1889 viennacl::linalg::ambm(lhs, lhs, T(1), 1,
false,
false, proxy.lhs(), proxy.rhs(), 1,
false,
false);
1894 template<
typename T,
typename ScalarType>
1895 struct op_executor<matrix_base<T>, op_inplace_sub, matrix_expression<const matrix_base<T>, const
ScalarType, op_mult> >
1897 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_base<T>,
const ScalarType, op_mult>
const & proxy)
1899 viennacl::linalg::ambm(lhs, lhs, T(1), 1,
false,
false, proxy.lhs(), proxy.rhs(), 1,
false,
true);
1907 template<
typename T,
typename LHS,
typename RHS,
typename OP,
typename ScalarType>
1908 struct op_executor<matrix_base<T>,
op_assign, matrix_expression<const matrix_expression<const LHS, const RHS, OP>, const
ScalarType, op_mult> >
1910 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<const LHS, const RHS, OP>,
const ScalarType, op_mult>
const & proxy)
1912 if (lhs.row_major())
1914 matrix<T> temp(proxy.lhs());
1915 lhs = temp * proxy.rhs();
1919 matrix<T, column_major> temp(proxy.lhs());
1920 lhs = temp * proxy.rhs();
1926 template<
typename T,
typename LHS,
typename RHS,
typename OP,
typename ScalarType>
1927 struct op_executor<matrix_base<T>, op_inplace_add, matrix_expression<const matrix_expression<const LHS, const RHS, OP>, const
ScalarType, op_mult> >
1929 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<const LHS, const RHS, OP>,
const ScalarType, op_mult>
const & proxy)
1931 if (lhs.row_major())
1933 matrix<T> temp(proxy.lhs());
1934 lhs += temp * proxy.rhs();
1938 matrix<T, column_major> temp(proxy.lhs());
1939 lhs += temp * proxy.rhs();
1945 template<
typename T,
typename LHS,
typename RHS,
typename OP,
typename ScalarType>
1946 struct op_executor<matrix_base<T>, op_inplace_sub, matrix_expression<const matrix_expression<const LHS, const RHS, OP>, const
ScalarType, op_mult> >
1948 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<const LHS, const RHS, OP>,
const ScalarType, op_mult>
const & proxy)
1950 if (lhs.row_major())
1952 matrix<T> temp(proxy.lhs());
1953 lhs -= temp * proxy.rhs();
1957 matrix<T, column_major> temp(proxy.lhs());
1958 lhs -= temp * proxy.rhs();
1967 template<
typename T,
typename ScalarType>
1968 struct op_executor<matrix_base<T>,
op_assign, matrix_expression<const matrix_base<T>, const
ScalarType, op_div> >
1970 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_base<T>,
const ScalarType, op_div>
const & proxy)
1977 template<
typename T,
typename ScalarType>
1978 struct op_executor<matrix_base<T>, op_inplace_add, matrix_expression<const matrix_base<T>, const
ScalarType, op_div> >
1980 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_base<T>,
const ScalarType, op_div>
const & proxy)
1982 viennacl::linalg::ambm(lhs, lhs, T(1), 1,
false,
false, proxy.lhs(), proxy.rhs(), 1,
true,
false);
1987 template<
typename T,
typename ScalarType>
1988 struct op_executor<matrix_base<T>, op_inplace_sub, matrix_expression<const matrix_base<T>, const
ScalarType, op_div> >
1990 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_base<T>,
const ScalarType, op_div>
const & proxy)
1992 viennacl::linalg::ambm(lhs, lhs, T(1), 1,
false,
false, proxy.lhs(), proxy.rhs(), 1,
true,
true);
2000 template<
typename T,
typename LHS,
typename RHS,
typename OP,
typename ScalarType>
2001 struct op_executor<matrix_base<T>,
op_assign, matrix_expression<const matrix_expression<const LHS, const RHS, OP>, const
ScalarType, op_div> >
2003 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<const LHS, const RHS, OP>,
const ScalarType, op_div>
const & proxy)
2005 if (lhs.row_major())
2007 matrix<T> temp(proxy.lhs());
2008 lhs = temp / proxy.rhs();
2012 matrix<T, column_major> temp(proxy.lhs());
2013 lhs = temp / proxy.rhs();
2019 template<
typename T,
typename LHS,
typename RHS,
typename OP,
typename ScalarType>
2020 struct op_executor<matrix_base<T>, op_inplace_add, matrix_expression<const matrix_expression<const LHS, const RHS, OP>, const
ScalarType, op_div> >
2022 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<const LHS, const RHS, OP>,
const ScalarType, op_div>
const & proxy)
2024 if (lhs.row_major())
2026 matrix<T> temp(proxy.lhs());
2027 lhs += temp / proxy.rhs();
2031 matrix<T, column_major> temp(proxy.lhs());
2032 lhs += temp / proxy.rhs();
2038 template<
typename T,
typename LHS,
typename RHS,
typename OP,
typename ScalarType>
2039 struct op_executor<matrix_base<T>, op_inplace_sub, matrix_expression<const matrix_expression<const LHS, const RHS, OP>, const
ScalarType, op_div> >
2041 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<const LHS, const RHS, OP>,
const ScalarType, op_div>
const & proxy)
2043 if (lhs.row_major())
2045 matrix<T, row_major> temp(proxy.lhs());
2046 lhs -= temp / proxy.rhs();
2050 matrix<T, column_major> temp(proxy.lhs());
2051 lhs -= temp / proxy.rhs();
2059 template<
typename T,
typename LHS,
typename RHS>
2060 struct op_executor<matrix_base<T>,
op_assign, matrix_expression<const LHS, const RHS, op_add> >
2063 template<
typename LHS1,
typename RHS1>
2064 static void apply(matrix_base<T> & lhs, matrix_expression<const LHS1, const RHS1, op_add>
const & proxy)
2066 bool op_aliasing_lhs =
op_aliasing(lhs, proxy.lhs());
2067 bool op_aliasing_rhs =
op_aliasing(lhs, proxy.rhs());
2069 if (op_aliasing_lhs || op_aliasing_rhs)
2071 matrix_base<T> temp(proxy.lhs());
2072 op_executor<matrix_base<T>, op_inplace_add, RHS>::apply(temp, proxy.rhs());
2077 op_executor<matrix_base<T>,
op_assign, LHS>::apply(lhs, proxy.lhs());
2078 op_executor<matrix_base<T>, op_inplace_add, RHS>::apply(lhs, proxy.rhs());
2083 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_base<T>,
const matrix_base<T>, op_add>
const & proxy)
2086 proxy.lhs(), T(1), 1,
false,
false,
2087 proxy.rhs(), T(1), 1,
false,
false);
2091 template<
typename ScalarType>
2092 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<
const matrix_base<T>,
const ScalarType, op_mult>,
2093 const matrix_base<T>,
2094 op_add>
const & proxy)
2097 proxy.lhs().lhs(), proxy.lhs().rhs(), 1,
false,
false,
2098 proxy.rhs(), T(1), 1,
false,
false);
2102 template<
typename ScalarType>
2103 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<
const matrix_base<T>,
const ScalarType, op_div>,
2104 const matrix_base<T>,
2105 op_add>
const & proxy)
2108 proxy.lhs().lhs(), proxy.lhs().rhs(), 1,
true,
false,
2109 proxy.rhs(), T(1), 1,
false,
false);
2113 template<
typename ScalarType>
2114 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_base<T>,
2115 const matrix_expression<
const matrix_base<T>,
const ScalarType, op_mult>,
2116 op_add>
const & proxy)
2119 proxy.lhs(), T(1), 1,
false,
false,
2120 proxy.rhs().lhs(), proxy.rhs().rhs(), 1,
false,
false);
2124 template<
typename ScalarType>
2125 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_base<T>,
2126 const matrix_expression<
const matrix_base<T>,
const ScalarType, op_div>,
2127 op_add>
const & proxy)
2130 proxy.lhs(), T(1), 1,
false,
false,
2131 proxy.rhs().lhs(), proxy.rhs().rhs(), 1,
true,
false);
2135 template<
typename ScalarType1,
typename ScalarType2>
2136 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<
const matrix_base<T>,
const ScalarType1, op_mult>,
2137 const matrix_expression<
const matrix_base<T>,
const ScalarType2, op_mult>,
2138 op_add>
const & proxy)
2141 proxy.lhs().lhs(), proxy.lhs().rhs(), 1,
false,
false,
2142 proxy.rhs().lhs(), proxy.rhs().rhs(), 1,
false,
false);
2146 template<
typename ScalarType1,
typename ScalarType2>
2147 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<
const matrix_base<T>,
const ScalarType1, op_mult>,
2148 const matrix_expression<
const matrix_base<T>,
const ScalarType2, op_div>,
2149 op_add>
const & proxy)
2152 proxy.lhs().lhs(), proxy.lhs().rhs(), 1,
false,
false,
2153 proxy.rhs().lhs(), proxy.rhs().rhs(), 1,
true,
false);
2157 template<
typename ScalarType1,
typename ScalarType2>
2158 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<
const matrix_base<T>,
const ScalarType1, op_div>,
2159 const matrix_expression<
const matrix_base<T>,
const ScalarType2, op_mult>,
2160 op_add>
const & proxy)
2163 proxy.lhs().lhs(), proxy.lhs().rhs(), 1,
true,
false,
2164 proxy.rhs().lhs(), proxy.rhs().rhs(), 1,
false,
false);
2168 template<
typename ScalarType1,
typename ScalarType2>
2169 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<
const matrix_base<T>,
const ScalarType1, op_div>,
2170 const matrix_expression<
const matrix_base<T>,
const ScalarType2, op_div>,
2171 op_add>
const & proxy)
2174 proxy.lhs().lhs(), proxy.lhs().rhs(), 1,
true,
false,
2175 proxy.rhs().lhs(), proxy.rhs().rhs(), 1,
true,
false);
2180 template<
typename T,
typename LHS,
typename RHS>
2181 struct op_executor<matrix_base<T>, op_assign, matrix_expression<const LHS, const RHS, op_prod> >
2183 template<
typename SparseMatrixType>
2184 static void apply(matrix_base<T> & lhs, matrix_expression<
const SparseMatrixType,
2191 matrix_base<T> temp(proxy);
2199 template<
typename SparseMatrixType >
2200 static void apply(matrix_base<T> & lhs, matrix_expression<
const SparseMatrixType,
2209 matrix_base<T> temp(proxy);
2219 template<
typename T,
typename LHS,
typename RHS>
2220 struct op_executor<matrix_base<T>, op_inplace_add, matrix_expression<const LHS, const RHS, op_add> >
2223 template<
typename LHS1,
typename RHS1>
2224 static void apply(matrix_base<T> & lhs, matrix_expression<const LHS1, const RHS1, op_add>
const & proxy)
2226 bool op_aliasing_lhs =
op_aliasing(lhs, proxy.lhs());
2227 bool op_aliasing_rhs =
op_aliasing(lhs, proxy.rhs());
2229 if (op_aliasing_lhs || op_aliasing_rhs)
2231 matrix_base<T> temp(proxy.lhs());
2232 op_executor<matrix_base<T>, op_inplace_add, RHS>::apply(temp, proxy.rhs());
2237 op_executor<matrix_base<T>, op_inplace_add, LHS>::apply(lhs, proxy.lhs());
2238 op_executor<matrix_base<T>, op_inplace_add, RHS>::apply(lhs, proxy.rhs());
2243 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_base<T>,
const matrix_base<T>, op_add>
const & proxy)
2246 proxy.lhs(), T(1), 1,
false,
false,
2247 proxy.rhs(), T(1), 1,
false,
false);
2251 template<
typename ScalarType>
2252 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<
const matrix_base<T>,
const ScalarType, op_mult>,
2253 const matrix_base<T>,
2254 op_add>
const & proxy)
2257 proxy.lhs().lhs(), proxy.lhs().rhs(), 1,
false,
false,
2258 proxy.rhs(), T(1), 1,
false,
false);
2262 template<
typename ScalarType>
2263 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<
const matrix_base<T>,
const ScalarType, op_div>,
2264 const matrix_base<T>,
2265 op_add>
const & proxy)
2268 proxy.lhs().lhs(), proxy.lhs().rhs(), 1,
true,
false,
2269 proxy.rhs(), T(1), 1,
false,
false);
2273 template<
typename ScalarType>
2274 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_base<T>,
2275 const matrix_expression<
const matrix_base<T>,
const ScalarType, op_mult>,
2276 op_add>
const & proxy)
2279 proxy.lhs(), T(1), 1,
false,
false,
2280 proxy.rhs().lhs(), proxy.rhs().rhs(), 1,
false,
false);
2284 template<
typename ScalarType>
2285 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_base<T>,
2286 const matrix_expression<
const matrix_base<T>,
const ScalarType, op_div>,
2287 op_add>
const & proxy)
2290 proxy.lhs(), T(1), 1,
false,
false,
2291 proxy.rhs().lhs(), proxy.rhs().rhs(), 1,
true,
false);
2295 template<
typename ScalarType1,
typename ScalarType2>
2296 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<
const matrix_base<T>,
const ScalarType1, op_mult>,
2297 const matrix_expression<
const matrix_base<T>,
const ScalarType2, op_mult>,
2298 op_add>
const & proxy)
2301 proxy.lhs().lhs(), proxy.lhs().rhs(), 1,
false,
false,
2302 proxy.rhs().lhs(), proxy.rhs().rhs(), 1,
false,
false);
2306 template<
typename ScalarType1,
typename ScalarType2>
2307 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<
const matrix_base<T>,
const ScalarType1, op_mult>,
2308 const matrix_expression<
const matrix_base<T>,
const ScalarType2, op_div>,
2309 op_add>
const & proxy)
2312 proxy.lhs().lhs(), proxy.lhs().rhs(), 1,
false,
false,
2313 proxy.rhs().lhs(), proxy.rhs().rhs(), 1,
true,
false);
2317 template<
typename ScalarType1,
typename ScalarType2>
2318 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<
const matrix_base<T>,
const ScalarType1, op_div>,
2319 const matrix_expression<
const matrix_base<T>,
const ScalarType2, op_mult>,
2320 op_add>
const & proxy)
2323 proxy.lhs().lhs(), proxy.lhs().rhs(), 1,
true,
false,
2324 proxy.rhs().lhs(), proxy.rhs().rhs(), 1,
false,
false);
2328 template<
typename ScalarType1,
typename ScalarType2>
2329 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<
const matrix_base<T>,
const ScalarType1, op_div>,
2330 const matrix_expression<
const matrix_base<T>,
const ScalarType2, op_div>,
2331 op_add>
const & proxy)
2334 proxy.lhs().lhs(), proxy.lhs().rhs(), 1,
true,
false,
2335 proxy.rhs().lhs(), proxy.rhs().rhs(), 1,
true,
false);
2342 template<
typename T,
typename LHS,
typename RHS>
2343 struct op_executor<matrix_base<T>, op_inplace_sub, matrix_expression<const LHS, const RHS, op_add> >
2346 template<
typename LHS1,
typename RHS1>
2347 static void apply(matrix_base<T> & lhs, matrix_expression<const LHS1, const RHS1, op_add>
const & proxy)
2349 bool op_aliasing_lhs =
op_aliasing(lhs, proxy.lhs());
2350 bool op_aliasing_rhs =
op_aliasing(lhs, proxy.rhs());
2352 if (op_aliasing_lhs || op_aliasing_rhs)
2354 matrix_base<T> temp(proxy.lhs());
2355 op_executor<matrix_base<T>, op_inplace_add, RHS>::apply(temp, proxy.rhs());
2360 op_executor<matrix_base<T>, op_inplace_sub, LHS>::apply(lhs, proxy.lhs());
2361 op_executor<matrix_base<T>, op_inplace_sub, RHS>::apply(lhs, proxy.rhs());
2366 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_base<T>,
const matrix_base<T>, op_add>
const & proxy)
2369 proxy.lhs(), T(1), 1,
false,
true,
2370 proxy.rhs(), T(1), 1,
false,
true);
2374 template<
typename ScalarType>
2375 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<
const matrix_base<T>,
const ScalarType, op_mult>,
2376 const matrix_base<T>,
2377 op_add>
const & proxy)
2380 proxy.lhs().lhs(), proxy.lhs().rhs(), 1,
false,
true,
2381 proxy.rhs(), T(1), 1,
false,
true);
2385 template<
typename ScalarType>
2386 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<
const matrix_base<T>,
const ScalarType, op_div>,
2387 const matrix_base<T>,
2388 op_add>
const & proxy)
2391 proxy.lhs().lhs(), proxy.lhs().rhs(), 1,
true,
true,
2392 proxy.rhs(), T(1), 1,
false,
true);
2396 template<
typename ScalarType>
2397 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_base<T>,
2398 const matrix_expression<
const matrix_base<T>,
const ScalarType, op_mult>,
2399 op_add>
const & proxy)
2402 proxy.lhs(), T(1), 1,
false,
true,
2403 proxy.rhs().lhs(), proxy.rhs().rhs(), 1,
false,
true);
2407 template<
typename ScalarType>
2408 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_base<T>,
2409 const matrix_expression<
const matrix_base<T>,
const ScalarType, op_div>,
2410 op_add>
const & proxy)
2413 proxy.lhs(), T(1), 1,
false,
true,
2414 proxy.rhs().lhs(), proxy.rhs().rhs(), 1,
true,
true);
2418 template<
typename ScalarType1,
typename ScalarType2>
2419 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<
const matrix_base<T>,
const ScalarType1, op_mult>,
2420 const matrix_expression<
const matrix_base<T>,
const ScalarType2, op_mult>,
2421 op_add>
const & proxy)
2424 proxy.lhs().lhs(), proxy.lhs().rhs(), 1,
false,
true,
2425 proxy.rhs().lhs(), proxy.rhs().rhs(), 1,
false,
true);
2429 template<
typename ScalarType1,
typename ScalarType2>
2430 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<
const matrix_base<T>,
const ScalarType1, op_mult>,
2431 const matrix_expression<
const matrix_base<T>,
const ScalarType2, op_div>,
2432 op_add>
const & proxy)
2435 proxy.lhs().lhs(), proxy.lhs().rhs(), 1,
false,
true,
2436 proxy.rhs().lhs(), proxy.rhs().rhs(), 1,
true,
true);
2440 template<
typename ScalarType1,
typename ScalarType2>
2441 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<
const matrix_base<T>,
const ScalarType1, op_div>,
2442 const matrix_expression<
const matrix_base<T>,
const ScalarType2, op_mult>,
2443 op_add>
const & proxy)
2446 proxy.lhs().lhs(), proxy.lhs().rhs(), 1,
true,
true,
2447 proxy.rhs().lhs(), proxy.rhs().rhs(), 1,
false,
true);
2451 template<
typename ScalarType1,
typename ScalarType2>
2452 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<
const matrix_base<T>,
const ScalarType1, op_div>,
2453 const matrix_expression<
const matrix_base<T>,
const ScalarType2, op_div>,
2454 op_add>
const & proxy)
2457 proxy.lhs().lhs(), proxy.lhs().rhs(), 1,
true,
true,
2458 proxy.rhs().lhs(), proxy.rhs().rhs(), 1,
true,
true);
2469 template<
typename T,
typename LHS,
typename RHS>
2470 struct op_executor<matrix_base<T>, op_assign, matrix_expression<const LHS, const RHS, op_sub> >
2473 template<
typename LHS1,
typename RHS1>
2474 static void apply(matrix_base<T> & lhs, matrix_expression<const LHS1, const RHS1, op_sub>
const & proxy)
2476 bool op_aliasing_lhs =
op_aliasing(lhs, proxy.lhs());
2477 bool op_aliasing_rhs =
op_aliasing(lhs, proxy.rhs());
2479 if (op_aliasing_lhs || op_aliasing_rhs)
2481 matrix_base<T> temp(proxy.lhs());
2482 op_executor<matrix_base<T>, op_inplace_sub, RHS>::apply(temp, proxy.rhs());
2487 op_executor<matrix_base<T>, op_assign, LHS>::apply(lhs, proxy.lhs());
2488 op_executor<matrix_base<T>, op_inplace_sub, RHS>::apply(lhs, proxy.rhs());
2493 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_base<T>,
const matrix_base<T>, op_sub>
const & proxy)
2496 proxy.lhs(), T(1), 1,
false,
false,
2497 proxy.rhs(), T(1), 1,
false,
true);
2501 template<
typename ScalarType>
2502 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<
const matrix_base<T>,
const ScalarType, op_mult>,
2503 const matrix_base<T>,
2504 op_sub>
const & proxy)
2507 proxy.lhs().lhs(), proxy.lhs().rhs(), 1,
false,
false,
2508 proxy.rhs(), T(1), 1,
false,
true);
2512 template<
typename ScalarType>
2513 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<
const matrix_base<T>,
const ScalarType, op_div>,
2514 const matrix_base<T>,
2515 op_sub>
const & proxy)
2518 proxy.lhs().lhs(), proxy.lhs().rhs(), 1,
true,
false,
2519 proxy.rhs(), T(1), 1,
false,
true);
2523 template<
typename ScalarType>
2524 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_base<T>,
2525 const matrix_expression<
const matrix_base<T>,
const ScalarType, op_mult>,
2526 op_sub>
const & proxy)
2529 proxy.lhs(), T(1), 1,
false,
false,
2530 proxy.rhs().lhs(), proxy.rhs().rhs(), 1,
false,
true);
2534 template<
typename ScalarType>
2535 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_base<T>,
2536 const matrix_expression<
const matrix_base<T>,
const ScalarType, op_div>,
2537 op_sub>
const & proxy)
2540 proxy.lhs(), T(1), 1,
false,
false,
2541 proxy.rhs().lhs(), proxy.rhs().rhs(), 1,
true,
true);
2545 template<
typename ScalarType1,
typename ScalarType2>
2546 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<
const matrix_base<T>,
const ScalarType1, op_mult>,
2547 const matrix_expression<
const matrix_base<T>,
const ScalarType2, op_mult>,
2548 op_sub>
const & proxy)
2551 proxy.lhs().lhs(), proxy.lhs().rhs(), 1,
false,
false,
2552 proxy.rhs().lhs(), proxy.rhs().rhs(), 1,
false,
true);
2556 template<
typename ScalarType1,
typename ScalarType2>
2557 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<
const matrix_base<T>,
const ScalarType1, op_mult>,
2558 const matrix_expression<
const matrix_base<T>,
const ScalarType2, op_div>,
2559 op_sub>
const & proxy)
2562 proxy.lhs().lhs(), proxy.lhs().rhs(), 1,
false,
false,
2563 proxy.rhs().lhs(), proxy.rhs().rhs(), 1,
true,
true);
2567 template<
typename ScalarType1,
typename ScalarType2>
2568 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<
const matrix_base<T>,
const ScalarType1, op_div>,
2569 const matrix_expression<
const matrix_base<T>,
const ScalarType2, op_mult>,
2570 op_sub>
const & proxy)
2573 proxy.lhs().lhs(), proxy.lhs().rhs(), 1,
true,
false,
2574 proxy.rhs().lhs(), proxy.rhs().rhs(), 1,
false,
true);
2578 template<
typename ScalarType1,
typename ScalarType2>
2579 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<
const matrix_base<T>,
const ScalarType1, op_div>,
2580 const matrix_expression<
const matrix_base<T>,
const ScalarType2, op_div>,
2581 op_sub>
const & proxy)
2584 proxy.lhs().lhs(), proxy.lhs().rhs(), 1,
true,
false,
2585 proxy.rhs().lhs(), proxy.rhs().rhs(), 1,
true,
true);
2591 template<
typename T,
typename LHS,
typename RHS>
2592 struct op_executor<matrix_base<T>, op_inplace_add, matrix_expression<const LHS, const RHS, op_sub> >
2595 template<
typename LHS1,
typename RHS1>
2596 static void apply(matrix_base<T> & lhs, matrix_expression<const LHS1, const RHS1, op_sub>
const & proxy)
2598 bool op_aliasing_lhs =
op_aliasing(lhs, proxy.lhs());
2599 bool op_aliasing_rhs =
op_aliasing(lhs, proxy.rhs());
2601 if (op_aliasing_lhs || op_aliasing_rhs)
2603 matrix_base<T> temp(proxy.lhs());
2604 op_executor<matrix_base<T>, op_inplace_sub, RHS>::apply(temp, proxy.rhs());
2609 op_executor<matrix_base<T>, op_inplace_add, LHS>::apply(lhs, proxy.lhs());
2610 op_executor<matrix_base<T>, op_inplace_sub, RHS>::apply(lhs, proxy.rhs());
2615 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_base<T>,
const matrix_base<T>, op_sub>
const & proxy)
2618 proxy.lhs(), T(1), 1,
false,
false,
2619 proxy.rhs(), T(1), 1,
false,
true);
2623 template<
typename ScalarType>
2624 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<
const matrix_base<T>,
const ScalarType, op_mult>,
2625 const matrix_base<T>,
2626 op_sub>
const & proxy)
2629 proxy.lhs().lhs(), proxy.lhs().rhs(), 1,
false,
false,
2630 proxy.rhs(), T(1), 1,
false,
true);
2634 template<
typename ScalarType>
2635 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<
const matrix_base<T>,
const ScalarType, op_div>,
2636 const matrix_base<T>,
2637 op_sub>
const & proxy)
2640 proxy.lhs().lhs(), proxy.lhs().rhs(), 1,
true,
false,
2641 proxy.rhs(), T(1), 1,
false,
true);
2645 template<
typename ScalarType>
2646 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_base<T>,
2647 const matrix_expression<
const matrix_base<T>,
const ScalarType, op_mult>,
2648 op_sub>
const & proxy)
2651 proxy.lhs(), T(1), 1,
false,
false,
2652 proxy.rhs().lhs(), proxy.rhs().rhs(), 1,
false,
true);
2656 template<
typename ScalarType>
2657 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_base<T>,
2658 const matrix_expression<
const matrix_base<T>,
const ScalarType, op_div>,
2659 op_sub>
const & proxy)
2662 proxy.lhs(), T(1), 1,
false,
false,
2663 proxy.rhs().lhs(), proxy.rhs().rhs(), 1,
true,
true);
2667 template<
typename ScalarType1,
typename ScalarType2>
2668 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<
const matrix_base<T>,
const ScalarType1, op_mult>,
2669 const matrix_expression<
const matrix_base<T>,
const ScalarType2, op_mult>,
2670 op_sub>
const & proxy)
2673 proxy.lhs().lhs(), proxy.lhs().rhs(), 1,
false,
false,
2674 proxy.rhs().lhs(), proxy.rhs().rhs(), 1,
false,
true);
2678 template<
typename ScalarType1,
typename ScalarType2>
2679 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<
const matrix_base<T>,
const ScalarType1, op_mult>,
2680 const matrix_expression<
const matrix_base<T>,
const ScalarType2, op_div>,
2681 op_sub>
const & proxy)
2684 proxy.lhs().lhs(), proxy.lhs().rhs(), 1,
false,
false,
2685 proxy.rhs().lhs(), proxy.rhs().rhs(), 1,
true,
true);
2689 template<
typename ScalarType1,
typename ScalarType2>
2690 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<
const matrix_base<T>,
const ScalarType1, op_div>,
2691 const matrix_expression<
const matrix_base<T>,
const ScalarType2, op_mult>,
2692 op_sub>
const & proxy)
2695 proxy.lhs().lhs(), proxy.lhs().rhs(), 1,
true,
false,
2696 proxy.rhs().lhs(), proxy.rhs().rhs(), 1,
false,
true);
2700 template<
typename ScalarType1,
typename ScalarType2>
2701 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<
const matrix_base<T>,
const ScalarType1, op_div>,
2702 const matrix_expression<
const matrix_base<T>,
const ScalarType2, op_div>,
2703 op_sub>
const & proxy)
2706 proxy.lhs().lhs(), proxy.lhs().rhs(), 1,
true,
false,
2707 proxy.rhs().lhs(), proxy.rhs().rhs(), 1,
true,
true);
2714 template<
typename T,
typename LHS,
typename RHS>
2715 struct op_executor<matrix_base<T>, op_inplace_sub, matrix_expression<const LHS, const RHS, op_sub> >
2718 template<
typename LHS1,
typename RHS1>
2719 static void apply(matrix_base<T> & lhs, matrix_expression<const LHS1, const RHS1, op_sub>
const & proxy)
2721 bool op_aliasing_lhs =
op_aliasing(lhs, proxy.lhs());
2722 bool op_aliasing_rhs =
op_aliasing(lhs, proxy.rhs());
2724 if (op_aliasing_lhs || op_aliasing_rhs)
2726 matrix_base<T> temp(proxy.lhs());
2727 op_executor<matrix_base<T>, op_inplace_sub, RHS>::apply(temp, proxy.rhs());
2732 op_executor<matrix_base<T>, op_inplace_sub, LHS>::apply(lhs, proxy.lhs());
2733 op_executor<matrix_base<T>, op_inplace_add, RHS>::apply(lhs, proxy.rhs());
2738 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_base<T>,
const matrix_base<T>, op_sub>
const & proxy)
2741 proxy.lhs(), T(1), 1,
false,
true,
2742 proxy.rhs(), T(1), 1,
false,
false);
2746 template<
typename ScalarType>
2747 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<
const matrix_base<T>,
const ScalarType, op_mult>,
2748 const matrix_base<T>,
2749 op_sub>
const & proxy)
2752 proxy.lhs().lhs(), proxy.lhs().rhs(), 1,
false,
true,
2753 proxy.rhs(), T(1), 1,
false,
false);
2757 template<
typename ScalarType>
2758 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<
const matrix_base<T>,
const ScalarType, op_div>,
2759 const matrix_base<T>,
2760 op_sub>
const & proxy)
2763 proxy.lhs().lhs(), proxy.lhs().rhs(), 1,
true,
true,
2764 proxy.rhs(), T(1), 1,
false,
false);
2768 template<
typename ScalarType>
2769 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_base<T>,
2770 const matrix_expression<
const matrix_base<T>,
const ScalarType, op_mult>,
2771 op_sub>
const & proxy)
2774 proxy.lhs(), T(1), 1,
false,
true,
2775 proxy.rhs().lhs(), proxy.rhs().rhs(), 1,
false,
false);
2779 template<
typename ScalarType>
2780 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_base<T>,
2781 const matrix_expression<
const matrix_base<T>,
const ScalarType, op_div>,
2782 op_sub>
const & proxy)
2785 proxy.lhs(), T(1), 1,
false,
true,
2786 proxy.rhs().lhs(), proxy.rhs().rhs(), 1,
true,
false);
2790 template<
typename ScalarType1,
typename ScalarType2>
2791 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<
const matrix_base<T>,
const ScalarType1, op_mult>,
2792 const matrix_expression<
const matrix_base<T>,
const ScalarType2, op_mult>,
2793 op_sub>
const & proxy)
2796 proxy.lhs().lhs(), proxy.lhs().rhs(), 1,
false,
true,
2797 proxy.rhs().lhs(), proxy.rhs().rhs(), 1,
false,
false);
2801 template<
typename ScalarType1,
typename ScalarType2>
2802 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<
const matrix_base<T>,
const ScalarType1, op_mult>,
2803 const matrix_expression<
const matrix_base<T>,
const ScalarType2, op_div>,
2804 op_sub>
const & proxy)
2807 proxy.lhs().lhs(), proxy.lhs().rhs(), 1,
false,
true,
2808 proxy.rhs().lhs(), proxy.rhs().rhs(), 1,
true,
false);
2812 template<
typename ScalarType1,
typename ScalarType2>
2813 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<
const matrix_base<T>,
const ScalarType1, op_div>,
2814 const matrix_expression<
const matrix_base<T>,
const ScalarType2, op_mult>,
2815 op_sub>
const & proxy)
2818 proxy.lhs().lhs(), proxy.lhs().rhs(), 1,
true,
true,
2819 proxy.rhs().lhs(), proxy.rhs().rhs(), 1,
false,
false);
2823 template<
typename ScalarType1,
typename ScalarType2>
2824 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<
const matrix_base<T>,
const ScalarType1, op_div>,
2825 const matrix_expression<
const matrix_base<T>,
const ScalarType2, op_div>,
2826 op_sub>
const & proxy)
2829 proxy.lhs().lhs(), proxy.lhs().rhs(), 1,
true,
true,
2830 proxy.rhs().lhs(), proxy.rhs().rhs(), 1,
true,
false);
2837 template<
typename T,
typename LHS>
2838 struct op_executor<matrix_base<T>, op_assign, matrix_expression<const LHS, const int, op_vector_diag> >
2840 static void apply(matrix_base<T> & lhs, matrix_expression<
const vector_base<T>,
const int, op_vector_diag>
const & proxy)
2847 template<
typename T,
typename LHS>
2848 struct op_executor<vector_base<T>, op_assign, vector_expression<const LHS, const int, op_matrix_diag> >
2850 static void apply(vector_base<T> & lhs, vector_expression<
const matrix_base<T>,
const int, op_matrix_diag>
const & proxy)
2856 template<
typename T,
typename LHS>
2857 struct op_executor<vector_base<T>, op_assign, vector_expression<const LHS, const unsigned int, op_row> >
2859 static void apply(vector_base<T> & lhs, vector_expression<
const matrix_base<T>,
const unsigned int, op_row>
const & proxy)
2866 template<
typename T,
typename LHS>
2867 struct op_executor<vector_base<T>, op_assign, vector_expression<const LHS, const unsigned int, op_column> >
2869 static void apply(vector_base<T> & lhs, vector_expression<
const matrix_base<T>,
const unsigned int, op_column>
const & proxy)
2877 template<
typename T>
2878 struct op_executor<vector_base<T>, op_assign, vector_expression<const matrix_base<T>, const matrix_base<T>, op_row_sum> >
2880 static void apply(vector_base<T> & lhs, vector_expression<
const matrix_base<T>,
const matrix_base<T>, op_row_sum>
const & proxy)
2886 template<
typename T,
typename LHS,
typename RHS,
typename OP>
2887 struct op_executor<vector_base<T>, op_assign, vector_expression<const matrix_expression<LHS, RHS, OP>, const matrix_expression<LHS, RHS, OP>, op_row_sum> >
2889 static void apply(vector_base<T> & lhs, vector_expression<
const matrix_expression<LHS, RHS, OP>,
const matrix_expression<LHS, RHS, OP>, op_row_sum>
const & proxy)
2891 matrix_base<T> tmp(proxy.lhs());
2896 template<
typename T>
2897 struct op_executor<vector_base<T>, op_assign, vector_expression<const matrix_base<T>, const matrix_base<T>, op_col_sum> >
2899 static void apply(vector_base<T> & lhs, vector_expression<
const matrix_base<T>,
const matrix_base<T>, op_col_sum>
const & proxy)
2906 template<
typename T,
typename LHS,
typename RHS,
typename OP>
2907 struct op_executor<vector_base<T>, op_assign, vector_expression<const matrix_expression<LHS, RHS, OP>, const matrix_expression<LHS, RHS, OP>, op_col_sum> >
2909 static void apply(vector_base<T> & lhs, vector_expression<
const matrix_expression<LHS, RHS, OP>,
const matrix_expression<LHS, RHS, OP>, op_col_sum>
const & proxy)
2911 matrix_base<T> tmp(proxy.lhs());
2919 template<
typename T,
typename LHS,
typename RHS,
typename OP>
2920 struct op_executor<matrix_base<T>, op_assign, matrix_expression<const LHS, const RHS, op_element_binary<OP> > >
2923 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_base<T>,
const matrix_base<T>, op_element_binary<OP> >
const & proxy)
2929 template<
typename LHS2,
typename RHS2,
typename OP2>
2930 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_base<T>,
const matrix_expression<const LHS2, const RHS2, OP2>, op_element_binary<OP> >
const & proxy)
2932 matrix_base<T> temp(proxy.rhs());
2937 template<
typename LHS1,
typename RHS1,
typename OP1>
2938 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<const LHS1, const RHS1, OP1>,
const matrix_base<T>, op_element_binary<OP> >
const & proxy)
2940 matrix_base<T> temp(proxy.lhs());
2945 template<
typename LHS1,
typename RHS1,
typename OP1,
2946 typename LHS2,
typename RHS2,
typename OP2>
2947 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<const LHS1, const RHS1, OP1>,
2948 const matrix_expression<const LHS2, const RHS2, OP2>,
2949 op_element_binary<OP> >
const & proxy)
2951 matrix_base<T> temp1(proxy.lhs());
2952 matrix_base<T> temp2(proxy.rhs());
2958 template<
typename T,
typename LHS,
typename RHS,
typename OP>
2959 struct op_executor<matrix_base<T>, op_inplace_add, matrix_expression<const LHS, const RHS, op_element_binary<OP> > >
2962 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_base<T>,
const matrix_base<T>, op_element_binary<OP> >
const & proxy)
2964 matrix_base<T> temp(proxy);
2969 template<
typename LHS2,
typename RHS2,
typename OP2>
2970 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_base<T>,
const matrix_expression<const LHS2, const RHS2, OP2>, op_element_binary<OP> >
const & proxy)
2972 matrix_base<T> temp(proxy.rhs());
2979 template<
typename LHS1,
typename RHS1,
typename OP1>
2980 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<const LHS1, const RHS1, OP1>,
const matrix_base<T>, op_element_binary<OP> >
const & proxy)
2982 matrix_base<T> temp(proxy.lhs());
2989 template<
typename LHS1,
typename RHS1,
typename OP1,
2990 typename LHS2,
typename RHS2,
typename OP2>
2991 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<const LHS1, const RHS1, OP1>,
2992 const matrix_expression<const LHS2, const RHS2, OP2>,
2993 op_element_binary<OP> >
const & proxy)
2995 matrix_base<T> temp1(proxy.lhs());
2996 matrix_base<T> temp2(proxy.rhs());
3004 template<
typename T,
typename LHS,
typename RHS,
typename OP>
3005 struct op_executor<matrix_base<T>, op_inplace_sub, matrix_expression<const LHS, const RHS, op_element_binary<OP> > >
3009 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_base<T>,
const matrix_base<T>, op_element_binary<OP> >
const & proxy)
3011 matrix_base<T> temp(proxy);
3016 template<
typename LHS2,
typename RHS2,
typename OP2>
3017 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_base<T>,
const matrix_expression<const LHS2, const RHS2, OP2>, op_element_binary<OP> >
const & proxy)
3019 matrix_base<T> temp(proxy.rhs());
3026 template<
typename LHS1,
typename RHS1,
typename OP1>
3027 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<const LHS1, const RHS1, OP1>,
const matrix_base<T>, op_element_binary<OP> >
const & proxy)
3029 matrix_base<T> temp(proxy.lhs());
3036 template<
typename LHS1,
typename RHS1,
typename OP1,
3037 typename LHS2,
typename RHS2,
typename OP2>
3038 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<const LHS1, const RHS1, OP1>,
3039 const matrix_expression<const LHS2, const RHS2, OP2>,
3040 op_element_binary<OP> >
const & proxy)
3042 matrix_base<T> temp1(proxy.lhs());
3043 matrix_base<T> temp2(proxy.rhs());
3052 template<
typename T,
typename LHS,
typename RHS,
typename OP>
3053 struct op_executor<matrix_base<T>, op_assign, matrix_expression<const LHS, const RHS, op_element_unary<OP> > >
3056 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_base<T>,
const matrix_base<T>, op_element_unary<OP> >
const & proxy)
3062 template<
typename LHS2,
typename RHS2,
typename OP2>
3063 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<const LHS2, const RHS2, OP2>,
3064 const matrix_expression<const LHS2, const RHS2, OP2>,
3065 op_element_unary<OP> >
const & proxy)
3067 matrix_base<T> temp(proxy.rhs());
3072 template<
typename T,
typename LHS,
typename RHS,
typename OP>
3073 struct op_executor<matrix_base<T>, op_inplace_add, matrix_expression<const LHS, const RHS, op_element_unary<OP> > >
3076 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_base<T>,
const matrix_base<T>, op_element_unary<OP> >
const & proxy)
3078 matrix_base<T> temp(proxy);
3083 template<
typename LHS2,
typename RHS2,
typename OP2>
3084 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<const LHS2, const RHS2, OP2>,
3085 const matrix_expression<const LHS2, const RHS2, OP2>,
3086 op_element_unary<OP> >
const & proxy)
3088 matrix_base<T> temp(proxy.rhs());
3094 template<
typename T,
typename LHS,
typename RHS,
typename OP>
3095 struct op_executor<matrix_base<T>, op_inplace_sub, matrix_expression<const LHS, const RHS, op_element_unary<OP> > >
3098 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_base<T>,
const matrix_base<T>, op_element_unary<OP> >
const & proxy)
3100 matrix_base<T> temp(proxy);
3105 template<
typename LHS2,
typename RHS2,
typename OP2>
3106 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<const LHS2, const RHS2, OP2>,
3107 const matrix_expression<const LHS2, const RHS2, OP2>,
3108 op_element_unary<OP> >
const & proxy)
3110 matrix_base<T> temp(proxy.rhs());
3121 template<
typename T>
3122 struct op_executor<matrix_base<T>, op_assign, matrix_expression<const matrix_base<T>, const matrix_base<T>, op_mat_mat_prod> >
3124 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_base<T>,
const matrix_base<T>, op_mat_mat_prod>
const & rhs)
3128 matrix_base<T> temp(rhs);
3137 template<
typename T>
3138 struct op_executor<matrix_base<T>, op_assign, matrix_expression<const matrix_base<T>,
3139 const matrix_expression<const matrix_base<T>, const matrix_base<T>, op_trans>,
3142 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_base<T>,
3143 const matrix_expression<
const matrix_base<T>,
const matrix_base<T>, op_trans>,
3144 op_mat_mat_prod>
const & rhs)
3148 matrix_base<T> temp(rhs);
3157 template<
typename T,
typename LhsT,
typename RhsT,
typename OpT>
3158 struct op_executor<matrix_base<T>,
3160 matrix_expression<const matrix_base<T>,
3161 const matrix_expression<const LhsT, const RhsT, OpT>,
3165 static void apply(matrix_base<T> & lhs,
3166 matrix_expression<
const matrix_base<T>,
3167 const matrix_expression<const LhsT, const RhsT, OpT>,
3168 op_mat_mat_prod>
const & rhs)
3170 matrix_base<T> temp(rhs.rhs());
3178 template<
typename T>
3179 struct op_executor<matrix_base<T>, op_assign, matrix_expression<const matrix_expression<const matrix_base<T>, const matrix_base<T>, op_trans>,
3180 const matrix_base<T>,
3183 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<
const matrix_base<T>,
const matrix_base<T>, op_trans>,
3184 const matrix_base<T>,
3185 op_mat_mat_prod>
const & rhs)
3189 matrix_base<T> temp(rhs);
3198 template<
typename T,
typename LhsT,
typename RhsT,
typename OpT>
3199 struct op_executor<matrix_base<T>,
3201 matrix_expression<const matrix_expression<const LhsT, const RhsT, OpT>,
3202 const matrix_base<T>,
3206 static void apply(matrix_base<T> & lhs,
3207 matrix_expression<
const matrix_expression<const LhsT, const RhsT, OpT>,
3208 const matrix_base<T>,
3209 op_mat_mat_prod>
const & rhs)
3211 matrix_base<T> temp(rhs.lhs());
3218 template<
typename T>
3219 struct op_executor<matrix_base<T>, op_assign, matrix_expression<const matrix_expression<const matrix_base<T>, const matrix_base<T>, op_trans>,
3220 const matrix_expression<const matrix_base<T>, const matrix_base<T>, op_trans>,
3223 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<
const matrix_base<T>,
const matrix_base<T>, op_trans>,
3224 const matrix_expression<
const matrix_base<T>,
const matrix_base<T>, op_trans>,
3225 op_mat_mat_prod>
const & rhs)
3229 matrix_base<T> temp(rhs);
3238 template<
typename T,
3239 typename LhsT1,
typename RhsT1,
typename OpT1,
3240 typename LhsT2,
typename RhsT2,
typename OpT2>
3241 struct op_executor<matrix_base<T>,
3243 matrix_expression<const matrix_expression<const LhsT1, const RhsT1, OpT1>,
3244 const matrix_expression<const LhsT2, const RhsT2, OpT2>,
3248 static void apply(matrix_base<T> & lhs,
3249 matrix_expression<
const matrix_expression<const LhsT1, const RhsT1, OpT1>,
3250 const matrix_expression<const LhsT2, const RhsT2, OpT2>,
3251 op_mat_mat_prod>
const & rhs)
3253 matrix_base<T> temp1(rhs.lhs());
3254 matrix_base<T> temp2(rhs.rhs());
3263 template<
typename T>
3264 struct op_executor<matrix_base<T>, op_inplace_add, matrix_expression<const matrix_base<T>, const matrix_base<T>, op_mat_mat_prod> >
3266 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_base<T>,
const matrix_base<T>, op_mat_mat_prod>
const & rhs)
3270 matrix_base<T> temp(rhs);
3279 template<
typename T>
3280 struct op_executor<matrix_base<T>, op_inplace_add, matrix_expression<const matrix_base<T>,
3281 const matrix_expression<const matrix_base<T>, const matrix_base<T>, op_trans>,
3284 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_base<T>,
3285 const matrix_expression<
const matrix_base<T>,
const matrix_base<T>, op_trans>,
3286 op_mat_mat_prod>
const & rhs)
3290 matrix_base<T> temp(rhs);
3299 template<
typename T,
typename LhsT,
typename RhsT,
typename OpT>
3300 struct op_executor<matrix_base<T>,
3302 matrix_expression<const matrix_base<T>,
3303 const matrix_expression<const LhsT, const RhsT, OpT>,
3307 static void apply(matrix_base<T> & lhs,
3308 matrix_expression<
const matrix_base<T>,
3309 const matrix_expression<const LhsT, const RhsT, OpT>,
3310 op_mat_mat_prod>
const & rhs)
3312 matrix_base<T> temp(rhs.rhs());
3319 template<
typename T>
3320 struct op_executor<matrix_base<T>, op_inplace_add, matrix_expression<const matrix_expression<const matrix_base<T>, const matrix_base<T>, op_trans>,
3321 const matrix_base<T>,
3324 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<
const matrix_base<T>,
const matrix_base<T>, op_trans>,
3325 const matrix_base<T>,
3326 op_mat_mat_prod>
const & rhs)
3330 matrix_base<T> temp(rhs);
3339 template<
typename T,
typename LhsT,
typename RhsT,
typename OpT>
3340 struct op_executor<matrix_base<T>,
3342 matrix_expression<const matrix_expression<const LhsT, const RhsT, OpT>,
3343 const matrix_base<T>,
3347 static void apply(matrix_base<T> & lhs,
3348 matrix_expression<
const matrix_expression<const LhsT, const RhsT, OpT>,
3349 const matrix_base<T>,
3350 op_mat_mat_prod>
const & rhs)
3352 matrix_base<T> temp(rhs.lhs());
3359 template<
typename T>
3360 struct op_executor<matrix_base<T>, op_inplace_add, matrix_expression<const matrix_expression<const matrix_base<T>, const matrix_base<T>, op_trans>,
3361 const matrix_expression<const matrix_base<T>, const matrix_base<T>, op_trans>,
3364 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<
const matrix_base<T>,
const matrix_base<T>, op_trans>,
3365 const matrix_expression<
const matrix_base<T>,
const matrix_base<T>, op_trans>,
3366 op_mat_mat_prod>
const & rhs)
3370 matrix_base<T> temp(rhs);
3380 template<
typename T,
3381 typename LhsT1,
typename RhsT1,
typename OpT1,
3382 typename LhsT2,
typename RhsT2,
typename OpT2>
3383 struct op_executor<matrix_base<T>,
3385 matrix_expression<const matrix_expression<const LhsT1, const RhsT1, OpT1>,
3386 const matrix_expression<const LhsT2, const RhsT2, OpT2>,
3390 static void apply(matrix_base<T> & lhs,
3391 matrix_expression<
const matrix_expression<const LhsT1, const RhsT1, OpT1>,
3392 const matrix_expression<const LhsT2, const RhsT2, OpT2>,
3393 op_mat_mat_prod>
const & rhs)
3395 matrix_base<T> temp1(rhs.lhs());
3396 matrix_base<T> temp2(rhs.rhs());
3404 template<
typename T>
3405 struct op_executor<matrix_base<T>, op_inplace_sub, matrix_expression<const matrix_base<T>, const matrix_base<T>, op_mat_mat_prod> >
3407 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_base<T>,
const matrix_base<T>, op_mat_mat_prod>
const & rhs)
3411 matrix_base<T> temp(rhs);
3420 template<
typename T>
3421 struct op_executor<matrix_base<T>, op_inplace_sub, matrix_expression<const matrix_base<T>,
3422 const matrix_expression<const matrix_base<T>, const matrix_base<T>, op_trans>,
3425 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_base<T>,
3426 const matrix_expression<
const matrix_base<T>,
const matrix_base<T>, op_trans>,
3427 op_mat_mat_prod>
const & rhs)
3431 matrix_base<T> temp(rhs);
3440 template<
typename T,
typename LhsT,
typename RhsT,
typename OpT>
3441 struct op_executor<matrix_base<T>,
3443 matrix_expression<const matrix_base<T>,
3444 const matrix_expression<const LhsT, const RhsT, OpT>,
3448 static void apply(matrix_base<T> & lhs,
3449 matrix_expression<
const matrix_base<T>,
3450 const matrix_expression<const LhsT, const RhsT, OpT>,
3451 op_mat_mat_prod>
const & rhs)
3453 matrix_base<T> temp(rhs.rhs());
3460 template<
typename T>
3461 struct op_executor<matrix_base<T>, op_inplace_sub, matrix_expression<const matrix_expression<const matrix_base<T>, const matrix_base<T>, op_trans>,
3462 const matrix_base<T>,
3465 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<
const matrix_base<T>,
const matrix_base<T>, op_trans>,
3466 const matrix_base<T>,
3467 op_mat_mat_prod>
const & rhs)
3471 matrix_base<T> temp(rhs);
3480 template<
typename T,
typename LhsT,
typename RhsT,
typename OpT>
3481 struct op_executor<matrix_base<T>,
3483 matrix_expression<const matrix_expression<const LhsT, const RhsT, OpT>,
3484 const matrix_base<T>,
3488 static void apply(matrix_base<T> & lhs,
3489 matrix_expression<
const matrix_expression<const LhsT, const RhsT, OpT>,
3490 const matrix_base<T>,
3491 op_mat_mat_prod>
const & rhs)
3493 matrix_base<T> temp(rhs.lhs());
3500 template<
typename T>
3501 struct op_executor<matrix_base<T>, op_inplace_sub, matrix_expression<const matrix_expression<const matrix_base<T>, const matrix_base<T>, op_trans>,
3502 const matrix_expression<const matrix_base<T>, const matrix_base<T>, op_trans>,
3505 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<
const matrix_base<T>,
const matrix_base<T>, op_trans>,
3506 const matrix_expression<
const matrix_base<T>,
const matrix_base<T>, op_trans>,
3507 op_mat_mat_prod>
const & rhs)
3511 matrix_base<T> temp(rhs);
3520 template<
typename T,
3521 typename LhsT1,
typename RhsT1,
typename OpT1,
3522 typename LhsT2,
typename RhsT2,
typename OpT2>
3523 struct op_executor<matrix_base<T>,
3525 matrix_expression<const matrix_expression<const LhsT1, const RhsT1, OpT1>,
3526 const matrix_expression<const LhsT2, const RhsT2, OpT2>,
3530 static void apply(matrix_base<T> & lhs,
3531 matrix_expression<
const matrix_expression<const LhsT1, const RhsT1, OpT1>,
3532 const matrix_expression<const LhsT2, const RhsT2, OpT2>,
3533 op_mat_mat_prod>
const & rhs)
3535 matrix_base<T> temp1(rhs.lhs());
3536 matrix_base<T> temp2(rhs.rhs());
3544 template<
typename T>
3545 struct op_executor<vector_base<T>, op_assign, vector_expression<const matrix_base<T>, const vector_base<T>, op_prod> >
3547 static void apply(vector_base<T> & lhs, vector_expression<
const matrix_base<T>,
const vector_base<T>, op_prod>
const & rhs)
3552 vector_base<T> temp(rhs);
3561 template<
typename T>
3562 struct op_executor<vector_base<T>, op_assign, vector_expression<const matrix_expression<const matrix_base<T>, const matrix_base<T>, op_trans>,
3563 const vector_base<T>,
3566 static void apply(vector_base<T> & lhs, vector_expression<
const matrix_expression<
const matrix_base<T>,
const matrix_base<T>, op_trans>,
3567 const vector_base<T>,
3568 op_prod>
const & rhs)
3573 vector_base<T> temp(rhs);
3582 template<
typename T,
typename LhsT,
typename RhsT,
typename OpT>
3583 struct op_executor<vector_base<T>,
3585 vector_expression<const matrix_expression<const LhsT, const RhsT, OpT>,
3586 const vector_base<T>,
3590 static void apply(vector_base<T> & lhs,
3591 vector_expression<
const matrix_expression<const LhsT, const RhsT, OpT>,
3592 const vector_base<T>,
3593 op_prod>
const & rhs)
3595 matrix_base<T> temp(rhs.lhs());
3601 template<
typename T,
typename LhsT,
typename RhsT,
typename OpT>
3602 struct op_executor<vector_base<T>,
3604 vector_expression<const matrix_base<T>,
3605 const vector_expression<const LhsT, const RhsT, OpT>,
3609 static void apply(vector_base<T> & lhs,
3610 vector_expression<
const matrix_base<T>,
3611 const vector_expression<const LhsT, const RhsT, OpT>,
3612 op_prod>
const & rhs)
3614 vector_base<T> x(rhs.rhs());
3620 template<
typename T,
3621 typename LhsT1,
typename RhsT1,
typename OpT1,
3622 typename LhsT2,
typename RhsT2,
typename OpT2>
3623 struct op_executor<vector_base<T>,
3625 vector_expression<const matrix_expression<const LhsT1, const RhsT1, OpT1>,
3626 const vector_expression<const LhsT2, const RhsT2, OpT2>,
3630 static void apply(vector_base<T> & lhs,
3631 vector_expression<
const matrix_expression<const LhsT1, const RhsT1, OpT1>,
3632 const vector_expression<const LhsT2, const RhsT2, OpT2>,
3633 op_prod>
const & rhs)
3635 matrix_base<T> A(rhs.lhs());
3636 vector_base<T> x(rhs.rhs());
3644 template<
typename T>
3645 struct op_executor<vector_base<T>, op_inplace_add, vector_expression<const matrix_base<T>, const vector_base<T>, op_prod> >
3647 static void apply(vector_base<T> & lhs, vector_expression<
const matrix_base<T>,
const vector_base<T>, op_prod>
const & rhs)
3649 vector_base<T> temp(rhs);
3655 template<
typename T>
3656 struct op_executor<vector_base<T>, op_inplace_add, vector_expression<const matrix_expression<const matrix_base<T>, const matrix_base<T>, op_trans>,
3657 const vector_base<T>,
3660 static void apply(vector_base<T> & lhs, vector_expression<
const matrix_expression<
const matrix_base<T>,
const matrix_base<T>, op_trans>,
3661 const vector_base<T>,
3662 op_prod>
const & rhs)
3664 vector_base<T> temp(rhs);
3670 template<
typename T,
typename LhsT,
typename RhsT,
typename OpT>
3671 struct op_executor<vector_base<T>,
3673 vector_expression<const matrix_expression<const LhsT, const RhsT, OpT>,
3674 const vector_base<T>,
3678 static void apply(vector_base<T> & lhs,
3679 vector_expression<
const matrix_expression<const LhsT, const RhsT, OpT>,
3680 const vector_base<T>,
3681 op_prod>
const & rhs)
3683 matrix_base<T> A(rhs.lhs());
3684 vector_base<T> y(lhs);
3691 template<
typename T,
typename LhsT,
typename RhsT,
typename OpT>
3692 struct op_executor<vector_base<T>,
3694 vector_expression<const matrix_base<T>,
3695 const vector_expression<const LhsT, const RhsT, OpT>,
3699 static void apply(vector_base<T> & lhs,
3700 vector_expression<
const matrix_base<T>,
3701 const vector_expression<const LhsT, const RhsT, OpT>,
3702 op_prod>
const & rhs)
3704 vector_base<T> x(rhs.rhs());
3705 vector_base<T> y(lhs);
3712 template<
typename T,
3713 typename LhsT1,
typename RhsT1,
typename OpT1,
3714 typename LhsT2,
typename RhsT2,
typename OpT2>
3715 struct op_executor<vector_base<T>,
3717 vector_expression<const matrix_expression<const LhsT1, const RhsT1, OpT1>,
3718 const vector_expression<const LhsT2, const RhsT2, OpT2>,
3722 static void apply(vector_base<T> & lhs,
3723 vector_expression<
const matrix_expression<const LhsT1, const RhsT1, OpT1>,
3724 const vector_expression<const LhsT2, const RhsT2, OpT2>,
3725 op_prod>
const & rhs)
3727 matrix_base<T> A(rhs.lhs());
3728 vector_base<T> x(rhs.rhs());
3729 vector_base<T> y(lhs);
3738 template<
typename T>
3739 struct op_executor<vector_base<T>, op_inplace_sub, vector_expression<const matrix_base<T>, const vector_base<T>, op_prod> >
3741 static void apply(vector_base<T> & lhs, vector_expression<
const matrix_base<T>,
const vector_base<T>, op_prod>
const & rhs)
3743 vector_base<T> temp(rhs);
3749 template<
typename T>
3750 struct op_executor<vector_base<T>, op_inplace_sub, vector_expression<const matrix_expression<const matrix_base<T>, const matrix_base<T>, op_trans>,
3751 const vector_base<T>,
3754 static void apply(vector_base<T> & lhs, vector_expression<
const matrix_expression<
const matrix_base<T>,
const matrix_base<T>, op_trans>,
3755 const vector_base<T>,
3756 op_prod>
const & rhs)
3758 vector_base<T> temp(rhs);
3764 template<
typename T,
typename LhsT,
typename RhsT,
typename OpT>
3765 struct op_executor<vector_base<T>,
3767 vector_expression<const matrix_expression<const LhsT, const RhsT, OpT>,
3768 const vector_base<T>,
3772 static void apply(vector_base<T> & lhs,
3773 vector_expression<
const matrix_expression<const LhsT, const RhsT, OpT>,
3774 const vector_base<T>,
3775 op_prod>
const & rhs)
3777 matrix_base<T> A(rhs.lhs());
3778 vector_base<T> y(lhs);
3785 template<
typename T,
typename LhsT,
typename RhsT,
typename OpT>
3786 struct op_executor<vector_base<T>,
3788 vector_expression<const matrix_base<T>,
3789 const vector_expression<const LhsT, const RhsT, OpT>,
3793 static void apply(vector_base<T> & lhs,
3794 vector_expression<
const matrix_base<T>,
3795 const vector_expression<const LhsT, const RhsT, OpT>,
3796 op_prod>
const & rhs)
3798 vector_base<T> x(rhs.rhs());
3799 vector_base<T> y(lhs);
3806 template<
typename T,
3807 typename LhsT1,
typename RhsT1,
typename OpT1,
3808 typename LhsT2,
typename RhsT2,
typename OpT2>
3809 struct op_executor<vector_base<T>,
3811 vector_expression<const matrix_expression<const LhsT1, const RhsT1, OpT1>,
3812 const vector_expression<const LhsT2, const RhsT2, OpT2>,
3816 static void apply(vector_base<T> & lhs,
3817 vector_expression<
const matrix_expression<const LhsT1, const RhsT1, OpT1>,
3818 const vector_expression<const LhsT2, const RhsT2, OpT2>,
3819 op_prod>
const & rhs)
3821 matrix_base<T> A(rhs.lhs());
3822 vector_base<T> x(rhs.rhs());
3823 vector_base<T> y(lhs);
3834 template<
typename T>
3835 struct op_executor<matrix_base<T>, op_assign, matrix_expression<const vector_base<T>, const vector_base<T>, op_prod> >
3837 static void apply(matrix_base<T> & lhs, matrix_expression<
const vector_base<T>,
const vector_base<T>, op_prod>
const & rhs)
3845 template<
typename T,
typename ScalarType>
3846 struct op_executor<matrix_base<T>, op_assign, matrix_expression< const matrix_expression<const vector_base<T>, const vector_base<T>, op_prod>,
3850 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<
const vector_base<T>,
const vector_base<T>, op_prod>,
3852 op_mult>
const & rhs)
3860 template<
typename T>
3861 struct op_executor<matrix_base<T>, op_inplace_add, matrix_expression<const vector_base<T>, const vector_base<T>, op_prod> >
3863 static void apply(matrix_base<T> & lhs, matrix_expression<
const vector_base<T>,
const vector_base<T>, op_prod>
const & rhs)
3870 template<
typename T,
typename ScalarType>
3871 struct op_executor<matrix_base<T>, op_inplace_add, matrix_expression< const matrix_expression<const vector_base<T>, const vector_base<T>, op_prod>,
3875 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<
const vector_base<T>,
const vector_base<T>, op_prod>,
3877 op_mult>
const & rhs)
3884 template<
typename T>
3885 struct op_executor<matrix_base<T>, op_inplace_sub, matrix_expression<const vector_base<T>, const vector_base<T>, op_prod> >
3887 static void apply(matrix_base<T> & lhs, matrix_expression<
const vector_base<T>,
const vector_base<T>, op_prod>
const & rhs)
3894 template<
typename T,
typename ScalarType>
3895 struct op_executor<matrix_base<T>, op_inplace_sub, matrix_expression< const matrix_expression<const vector_base<T>, const vector_base<T>, op_prod>,
3899 static void apply(matrix_base<T> & lhs, matrix_expression<
const matrix_expression<
const vector_base<T>,
const vector_base<T>, op_prod>,
3901 op_mult>
const & rhs)
void row_sum_impl(const matrix_base< NumericT > &A, vector_base< NumericT > &result)
Simple enable-if variant that uses the SFINAE pattern.
A tag class representing multiplication by a scalar.
Helper class implementing an array on the host. Default case: No conversion necessary.
void matrix_diag_to_vector(const matrix_base< NumericT > &A, int k, vector_base< NumericT > &v)
Dispatcher interface for v = diag(A, k)
static vcl_size_t mem_index(vcl_size_t i, vcl_size_t j, vcl_size_t, vcl_size_t num_cols)
Returns the memory offset for entry (i,j) of a dense matrix.
void memory_write(mem_handle &dst_buffer, vcl_size_t dst_offset, vcl_size_t bytes_to_write, const void *ptr, bool async=false)
Writes data from main RAM identified by 'ptr' to the buffer identified by 'dst_buffer'.
viennacl::enable_if< viennacl::is_any_sparse_matrix< M1 >::value, matrix_expression< const M1, const M1, op_trans > >::type trans(const M1 &mat)
Returns an expression template class representing a transposed matrix.
Worker class for decomposing expression templates.
A proxy class for a single element of a vector or matrix. This proxy should not be noticed by end-use...
Implementations of dense matrix related operations including matrix-vector products.
size_type internal_size() const
Returns the total amount of allocated memory in multiples of sizeof(NumericT)
Class for representing strided submatrices of a bigger matrix A.
self_type & operator=(const self_type &other)
Helper class for checking whether a matrix has a row-major layout.
Helper struct for checking whether a type represents a sign flip on a viennacl::scalar<> ...
void matrix_assign(matrix_base< NumericT > &mat, NumericT s, bool clear=false)
matrix_expression< const self_type, const NumericT, op_mult > operator-() const
Sign flip for the matrix. Emulated to be equivalent to -1.0 * matrix.
A tag class representing the extraction of a matrix column to a vector.
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.)
A tag class representing a matrix given by a vector placed on a certain (off-)diagonal.
A tag class representing subtraction.
void am(matrix_base< NumericT > &mat1, matrix_base< NumericT > const &mat2, ScalarType1 const &alpha, vcl_size_t len_alpha, bool reciprocal_alpha, bool flip_sign_alpha)
viennacl::context context() const
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...
A tag indicating iteration along increasing row index of a matrix.
Expression template class for representing a tree of expressions which ultimately result in a matrix...
void ambm(matrix_base< NumericT > &mat1, matrix_base< NumericT > const &mat2, ScalarType1 const &alpha, vcl_size_t len_alpha, bool reciprocal_alpha, bool flip_sign_alpha, matrix_base< NumericT > const &mat3, ScalarType2 const &beta, vcl_size_t len_beta, bool reciprocal_beta, bool flip_sign_beta)
void clear(VectorType &vec)
Generic routine for setting all entries of a vector to zero. This is the version for non-ViennaCL obj...
This file provides the forward declarations for the main types used within ViennaCL.
base_type & operator=(viennacl::matrix_slice< viennacl::matrix< OtherNumericT, F2 > > const &B)
A tag class representing division.
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.
entry_proxy< NumericT > operator()(size_type row_index, size_type col_index)
Read-write access to a single element of the matrix/matrix_range/matrix_slice.
matrix_iterator(MatrixT &mat, vcl_size_t start_row, vcl_size_t start_col)
Represents a vector consisting of 1 at a given index and zeros otherwise. To be used as an initialize...
viennacl::enable_if< viennacl::is_scalar< S1 >::value, matrix_base< NumericT > & >::type operator/=(matrix_base< NumericT > &m1, S1 const &gpu_val)
Scales a matrix by a GPU scalar value.
viennacl::enable_if< viennacl::is_any_scalar< S1 >::value, matrix_expression< const matrix_base< NumericT >, const S1, op_mult >>::type operator*(S1 const &value, matrix_base< NumericT > const &m1)
Operator overload for the expression alpha * m1, where alpha is a host scalar (float or double) and m...
viennacl::scalar< float > s1
vcl_size_t internal_size(vector_base< NumericT > const &vec)
Helper routine for obtaining the buffer length of a ViennaCL vector.
An expression template class that represents a binary operation that yields a vector.
void element_op(matrix_base< T > &A, matrix_expression< const matrix_base< T >, const matrix_base< T >, OP > const &proxy)
Implementation of the element-wise operation A = B .* C and A = B ./ C for matrices (using MATLAB syn...
Forward declaration of dense matrix classes.
bool op_aliasing(vector_base< NumericT > const &, B const &)
result_of::size_type< MatrixType >::type size2(MatrixType const &mat)
Generic routine for obtaining the number of columns of a matrix (ViennaCL, uBLAS, etc...
matrix(size_type rows, size_type columns, viennacl::context ctx=viennacl::context())
Creates the matrix with the given dimensions.
vcl_size_t size1() const
Returns the size of the result vector.
viennacl::vector< NumericT > operator-(const vector_base< NumericT > &v1, const vector_expression< const matrix_base< NumericT >, const vector_base< NumericT >, op_prod > &proxy)
Implementation of the operation 'result = v1 - A * v2', where A is a matrix.
self_type & operator++(void)
self_type operator++(int)
A tag class representing the (off-)diagonal of a matrix.
base_type::size_type size_type
Represents a generic 'context' similar to an OpenCL context, but is backend-agnostic and thus also su...
Main namespace in ViennaCL. Holds all the basic types such as vector, matrix, etc. and defines operations upon them.
value_type operator*(void)
viennacl::enable_if< viennacl::is_scalar< S1 >::value, matrix_base< NumericT > & >::type operator*=(matrix_base< NumericT > &m1, S1 const &gpu_val)
Scales a matrix by a GPU scalar value.
matrix(zero_matrix< NumericT > const &m)
Creates the matrix from the supplied zero matrix.
vcl_size_t size(VectorType const &vec)
Generic routine for obtaining the size of a vector (ViennaCL, uBLAS, etc.)
matrix(NumericT *ptr_to_mem, viennacl::memory_types mem_type, size_type rows, size_type internal_row_count, size_type cols, size_type internal_col_count)
Wraps a CUDA or host buffer provided by the user including padding of rows and columns.
matrix_base()
The default constructor. Does not allocate any memory.
void convert(matrix_base< DestNumericT > &dest, matrix_base< SrcNumericT > const &src)
void matrix_row(const matrix_base< NumericT > &A, unsigned int i, vector_base< NumericT > &v)
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.
void clear()
Resets all entries to zero.
Implementations of operations using sparse matrices.
A tag class representing addition.
Represents a vector consisting of scalars 's' only, i.e. v[i] = s for all i. To be used as an initial...
matrix(const self_type &other)
matrix_expression(LHS &lhs, RHS &rhs)
void scaled_rank_1_update(matrix_base< NumericT > &mat1, S1 const &alpha, vcl_size_t len_alpha, bool reciprocal_alpha, bool flip_sign_alpha, const vector_base< NumericT > &vec1, const vector_base< NumericT > &vec2)
The implementation of the operation mat += alpha * vec1 * vec2^T, i.e. a scaled rank 1 update...
viennacl::enable_if< viennacl::is_any_scalar< S1 >::value, matrix_expression< const matrix_expression< const LHS, const RHS, OP >, const S1, op_div > >::type operator/(matrix_expression< const LHS, const RHS, OP > const &proxy, S1 const &val)
Operator overload for the division of a matrix expression by a scalar from the right, e.g. (beta * m1) / alpha. Here, beta * m1 is wrapped into a matrix_expression and then divided by alpha.
Determines whether a given expression has a row-major matrix layout.
size_type size2() const
Returns the number of columns.
handle_type & handle()
Returns the OpenCL handle, non-const-version.
matrix(const base_type &other)
void trans(const matrix_expression< const matrix_base< NumericT, SizeT, DistanceT >, const matrix_base< NumericT, SizeT, DistanceT >, op_trans > &proxy, matrix_base< NumericT > &temp_trans)
vector_expression< const matrix_base< NumericT >, const int, op_matrix_diag > diag(const matrix_base< NumericT > &A, int k=0)
viennacl::memory_types active_handle_id(T const &obj)
Returns an ID for the currently active memory domain of an object.
Represents a vector consisting of zeros only. To be used as an initializer for viennacl::vector, vector_range, or vector_slize only.
void matrix_column(const matrix_base< NumericT > &A, unsigned int j, vector_base< NumericT > &v)
base_type & operator=(viennacl::matrix_range< viennacl::matrix< OtherNumericT, F2 > > const &B)
size_type size1() const
Returns the number of rows.
vector_expression< const matrix_base< NumericT, F >, const unsigned int, op_row > row(const matrix_base< NumericT, F > &A, unsigned int i)
RHS & rhs() const
Get right hand side operand.
viennacl::vector< NumericT > operator+(const vector_base< NumericT > &v1, const vector_expression< const matrix_base< NumericT >, const vector_base< NumericT >, op_prod > &proxy)
Implementation of the operation 'result = v1 + A * v2', where A is a matrix.
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...
A tag class representing matrix-vector products and element-wise multiplications. ...
viennacl::context context(T const &t)
Returns an ID for the currently active memory domain of an object.
matrix(matrix_expression< LHS, RHS, OP > const &proxy)
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) ...
bool row_major(T const &)
void set(vcl_size_t index, U value)
matrix(scalar_matrix< NumericT > const &m)
Creates the matrix from the supplied scalar matrix.
void resize(size_type rows, size_type columns, bool preserve=true)
Resizes the matrix. Existing entries can optionally be preserved.
A tag class representing transposed matrices.
Helper implementations that deduce the dimensions of the supplied matrix-valued expressions.
vcl_size_t raw_size() const
Returns the number of bytes of the currently active buffer.
size_type internal_size2() const
Returns the internal number of columns. Usually required for launching OpenCL kernels only...
Class for representing non-strided submatrices of a bigger matrix A.
self_type & operator*=(char val)
Scales the matrix by a char (8-bit integer)
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...
viennacl::matrix< float > m1
MatrixT::value_type value_type
void matrix_diagonal_assign(matrix_base< NumericT > &mat, NumericT s)
size_type internal_size1() const
Returns the internal number of rows. Usually required for launching OpenCL kernels only...
bool operator!=(self_type const &other)
static vcl_size_t mem_index(vcl_size_t i, vcl_size_t j, vcl_size_t num_rows, vcl_size_t)
Returns the memory offset for entry (i,j) of a dense matrix.
vector_expression< const matrix_base< NumericT, F >, const unsigned int, op_column > column(const matrix_base< NumericT, F > &A, unsigned int j)
Extracts the underlying OpenCL handle from a vector, a matrix, an expression etc. ...
self_type & operator+=(const matrix_expression< const LHS, const RHS, OP > &proxy)
matrix(NumericT *ptr_to_mem, viennacl::memory_types mem_type, size_type rows, size_type cols)
Wraps a CUDA or host buffer provided by the user.
self_type & operator/=(char val)
Scales the matrix by a char (8-bit integer)
void prod_impl(const matrix_base< NumericT > &mat, const vector_base< NumericT > &vec, vector_base< NumericT > &result)
Carries out matrix-vector multiplication.
void resize(size_type rows, size_type columns, bool preserve=true)
viennacl::backend::mem_handle & handle(T &obj)
Returns the generic memory handle of an object. Non-const version.
matrix(identity_matrix< NumericT > const &m)
Creates the matrix from the supplied identity matrix.
void matrix_diag_from_vector(const vector_base< NumericT > &v, int k, matrix_base< NumericT > &A)
Dispatcher interface for A = diag(v, k)
LHS & lhs() const
Get left hand side operand.
A tag class representing the extraction of a matrix row to a vector.
A proxy class for a single element of a vector or matrix. This proxy should not be noticed by end-use...
MatrixT & operator()(void) const
void ambm_m(matrix_base< NumericT > &mat1, matrix_base< NumericT > const &mat2, ScalarType1 const &alpha, vcl_size_t len_alpha, bool reciprocal_alpha, bool flip_sign_alpha, matrix_base< NumericT > const &mat3, ScalarType2 const &beta, vcl_size_t len_beta, bool reciprocal_beta, bool flip_sign_beta)
Implementation of the ViennaCL scalar class.
A collection of compile time type deductions.
A tag for row-major storage of a dense matrix.
matrix()
The default constructor. Does not allocate any memory.
bool operator==(self_type const &other)
ram_handle_type & ram_handle()
Returns the handle to a buffer in CPU RAM. NULL is returned if no such buffer has been allocated...
A tag indicating iteration along increasing columns index of a matrix.
Simple enable-if variant that uses the SFINAE pattern.
void column_sum_impl(const matrix_base< NumericT > &A, vector_base< NumericT > &result)
self_type & operator-=(const matrix_expression< const LHS, const RHS, OP > &proxy)
base_type & operator=(viennacl::matrix< OtherNumericT, F2 > const &B)
void fast_copy(const const_vector_iterator< SCALARTYPE, ALIGNMENT > &gpu_begin, const const_vector_iterator< SCALARTYPE, ALIGNMENT > &gpu_end, CPU_ITERATOR cpu_begin)