1 #ifndef _theplu_yat_utility_blas_level2 2 #define _theplu_yat_utility_blas_level2 25 #include "BasicMatrix.h" 26 #include "BasicVector.h" 27 #include "BLAS_level3.h" 28 #include "VectorExpression.h" 29 #include "yat_assert.h" 31 #include "expression/MatrixTraits.h" 33 #include "gsl/gsl_blas.h" 46 namespace expression {
48 template<
typename MATRIX,
typename VECTOR>
49 class MatrixVector :
public VectorExpression<MatrixVector<MATRIX, VECTOR> >
53 MatrixVector(
const BasicMatrix<MATRIX>& lhs,
54 const BasicVector<VECTOR>& rhs)
56 YAT_ASSERT(lhs.columns() == rhs.size());
57 YAT_ASSERT(lhs.rows());
58 this->allocate_memory(lhs.rows());
60 MatrixTraits<BasicMatrix<MATRIX> > traits;
62 gsl_blas_dgemv(traits.transpose_type(), traits.factor(lhs),
63 traits.get_gsl_matrix_p(lhs),
64 rhs.gsl_vector_p(), 0.0, this->v_);
69 double operator()(
size_t i)
const 72 return gsl_vector_get(this->v_, i);
76 size_t size(
void)
const 79 return this->v_->size;
83 void calculate_gsl_vector_p(
void)
const 105 template<
class MATRIX,
class VECTOR>
106 expression::MatrixVector<MATRIX, VECTOR>
110 return expression::MatrixVector<MATRIX, VECTOR>(lhs, rhs);
122 template<
class MATRIX,
class VECTOR>
124 expression::TransposedMatrix<MATRIX> >,
128 YAT_ASSERT(lhs.
size() == rhs.
rows());
129 return transpose(rhs) * lhs;
size_t size(void) const
Definition: BasicVector.h:71
expression::MatrixVector< MatrixExpression< expression::TransposedMatrix< MATRIX > >, VECTOR > operator*(const BasicVector< VECTOR > &lhs, const BasicMatrix< MATRIX > &rhs)
Definition: BLAS_level2.h:126
The Department of Theoretical Physics namespace as we define it.
size_t columns(void) const
Definition: BasicMatrix.h:67
size_t rows(void) const
Definition: BasicMatrix.h:61
An expression that can be converted to a Matrix.
Definition: MatrixExpression.h:46
Definition: BasicVector.h:48
Definition: BasicMatrix.h:38
expression::MatrixVector< MATRIX, VECTOR > operator*(const BasicMatrix< MATRIX > &lhs, const BasicVector< VECTOR > &rhs)
Definition: BLAS_level2.h:107