Wiki

Clone wiki

blaze / BLAS Functions


For vector/vector, matrix/vector and matrix/matrix multiplications with large dense matrices Blaze relies on the efficiency of BLAS libraries. For this purpose, Blaze implements several convenient C++ wrapper functions for several BLAS functions. The following sections give a complete overview of all available BLAS level 1, 2 and 3 functions.


BLAS Level 1

Dot Product (dotu)

The following wrapper functions provide a generic interface for the BLAS functions for the dot product of two dense vectors (cblas_sdot(), cblas_ddot(), cblas_cdotu_sub(), and cblas_zdotu_sub()):

namespace blaze {

float dotu( blas_int_t n, const float* x, blas_int_t incX, const float* y, blas_int_t incY );

double dotu( blas_int_t n, const double* x, blas_int_t incX, const double* y, blas_int_t incY );

complex<float> dotu( blas_int_t n, const complex<float>* x, blas_int_t incX,
                     const complex<float>* y, blas_int_t incY );

complex<double> dotu( blas_int_t n, const complex<double>* x, blas_int_t incX,
                      const complex<double>* y, blas_int_t incY );

template< typename VT1, bool TF1, typename VT2, bool TF2 >
ElementType_<VT1> dotu( const DenseVector<VT1,TF1>& x, const DenseVector<VT2,TF2>& y );

} // namespace blaze

Complex Conjugate Dot Product (dotc)

The following wrapper functions provide a generic interface for the BLAS functions for the complex conjugate dot product of two dense vectors (cblas_sdot(), cblas_ddot(), cblas_cdotc_sub(), and cblas_zdotc_sub()):

namespace blaze {

float dotc( blas_int_t n, const float* x, blas_int_t incX, const float* y, blas_int_t incY );

double dotc( blas_int_t n, const double* x, blas_int_t incX, const double* y, blas_int_t incY );

complex<float> dotc( blas_int_t n, const complex<float>* x, blas_int_t incX,
                     const complex<float>* y, blas_int_t incY );

complex<double> dotc( blas_int_t n, const complex<double>* x, blas_int_t incX,
                      const complex<double>* y, blas_int_t incY );

template< typename VT1, bool TF1, typename VT2, bool TF2 >
ElementType_<VT1> dotc( const DenseVector<VT1,TF1>& x, const DenseVector<VT2,TF2>& y );

} // namespace blaze

Axpy Product (axpy)

The following wrapper functions provide a generic interface for the BLAS functions for the axpy product of two dense vectors (cblas_saxpy(), cblas_daxpy(), cblas_caxpy(), and cblas_zaxpy()):

namespace blaze {

void axpy( blas_int_t n, float alpha, const float* x, blas_int_t incX, float* y, blas_int_t incY );

void axpy( blas_int_t n, double alpha, const double* x, blas_int_t incX, double* y, blas_int_t incY );

void axpy( blas_int_t n, complex<float> alpha, const complex<float>* x,
           blas_int_t incX, complex<float>* y, blas_int_t incY );

void axpy( blas_int_t n, complex<double> alpha, const complex<double>* x,
           blas_int_t incX, complex<double>* y, blas_int_t incY );

template< typename VT1, bool TF1, typename VT2, bool TF2, typename ST >
void axpy( const DenseVector<VT1,TF1>& x, const DenseVector<VT2,TF2>& y, ST alpha );

} // namespace blaze

BLAS Level 2

General Matrix/Vector Multiplication (gemv)

The following wrapper functions provide a generic interface for the BLAS functions for the general matrix/vector multiplication (cblas_sgemv(), cblas_dgemv(), cblas_cgemv(), and cblas_zgemv()):

namespace blaze {

void gemv( CBLAS_ORDER layout, CBLAS_TRANSPOSE transA, blas_int_t m, blas_int_t n,
           float alpha, const float* A, blas_int_t lda, const float* x, blas_int_t incX,
           float beta, float* y, blas_int_t incY );

void gemv( CBLAS_ORDER layout, CBLAS_TRANSPOSE transA, blas_int_t m, blas_int_t n,
           double alpha, const double* A, blas_int_t lda, const double* x, blas_int_t incX,
           double beta, double* y, blas_int_t incY );

void gemv( CBLAS_ORDER layout, CBLAS_TRANSPOSE transA, blas_int_t m, blas_int_t n,
           complex<float> alpha, const complex<float>* A, blas_int_t lda,
           const complex<float>* x, blas_int_t incX, complex<float> beta,
           complex<float>* y, blas_int_t incY );

void gemv( CBLAS_ORDER layout, CBLAS_TRANSPOSE transA, blas_int_t m, blas_int_t n,
           complex<double> alpha, const complex<double>* A, blas_int_t lda,
           const complex<double>* x, blas_int_t incX, complex<double> beta,
           complex<double>* y, blas_int_t incY );

} // namespace blaze

Triangular Matrix/Vector Multiplication (trmv)

The following wrapper functions provide a generic interface for the BLAS functions for the matrix/vector multiplication with a triangular matrix (cblas_strmv(), cblas_dtrmv(), cblas_ctrmv(), and cblas_ztrmv()):

namespace blaze {

void trmv( CBLAS_ORDER order, CBLAS_UPLO uplo, CBLAS_TRANSPOSE transA, CBLAS_DIAG diag,
           blas_int_t n, const float* A, blas_int_t lda, float* x, blas_int_t incX );

void trmv( CBLAS_ORDER order, CBLAS_UPLO uplo, CBLAS_TRANSPOSE transA, CBLAS_DIAG diag,
           blas_int_t n, const double* A, blas_int_t lda, double* x, blas_int_t incX );

void trmv( CBLAS_ORDER order, CBLAS_UPLO uplo, CBLAS_TRANSPOSE transA, CBLAS_DIAG diag,
           blas_int_t n, const complex<float>* A, blas_int_t lda, complex<float>* x, blas_int_t incX );

void trmv( CBLAS_ORDER order, CBLAS_UPLO uplo, CBLAS_TRANSPOSE transA, CBLAS_DIAG diag,
           blas_int_t n, const complex<double>* A, blas_int_t lda, complex<double>* x, blas_int_t incX );

template< typename VT, typename MT, bool SO >
void trmv( DenseVector<VT,false>& x, const DenseMatrix<MT,SO>& A, CBLAS_UPLO uplo );

template< typename VT, typename MT, bool SO >
void trmv( DenseVector<VT,true>& x, const DenseMatrix<MT,SO>& A, CBLAS_UPLO uplo );

} // namespace blaze

BLAS Level 3

General Matrix/Matrix Multiplication (gemm)

The following wrapper functions provide a generic interface for the BLAS functions for the general matrix/matrix multiplication (cblas_sgemm(), cblas_dgemm(), cblas_cgemm(), and cblas_zgemm()):

namespace blaze {

void gemm( CBLAS_ORDER order, CBLAS_TRANSPOSE transA, CBLAS_TRANSPOSE transB,
           blas_int_t m, blas_int_t n, blas_int_t k, float alpha, const float* A,
           blas_int_t lda, const float* B, blas_int_t ldb, float beta, float* C,
           blas_int_t ldc );

void gemm( CBLAS_ORDER order, CBLAS_TRANSPOSE transA, CBLAS_TRANSPOSE transB,
           blas_int_t m, blas_int_t n, blas_int_t k, double alpha, const double* A,
           blas_int_t lda, const double* B, blas_int_t ldb, double beta, float* C,
           blas_int_t ldc );

void gemm( CBLAS_ORDER order, CBLAS_TRANSPOSE transA, CBLAS_TRANSPOSE transB,
           blas_int_t m, blas_int_t n, blas_int_t k, complex<float> alpha,
           const complex<float>* A, blas_int_t lda, const complex<float>* B,
           blas_int_t ldb, complex<float> beta, float* C, blas_int_t ldc );

void gemm( CBLAS_ORDER order, CBLAS_TRANSPOSE transA, CBLAS_TRANSPOSE transB,
           blas_int_t m, blas_int_t n, blas_int_t k, complex<double> alpha,
           const complex<double>* A, blas_int_t lda, const complex<double>* B,
           blas_int_t ldb, complex<double> beta, float* C, blas_int_t ldc );x

} // namespace blaze

Triangular Matrix/Matrix Multiplication (trmm)

The following wrapper functions provide a generic interface for the BLAS functions for the matrix/matrix multiplication with a triangular matrix (cblas_strmm(), cblas_dtrmm(), cblas_ctrmm(), and cblas_ztrmm()):

namespace blaze {

void trmm( CBLAS_ORDER order, CBLAS_SIDE side, CBLAS_UPLO uplo, CBLAS_TRANSPOSE transA,
           CBLAS_DIAG diag, blas_int_t m, blas_int_t n, float alpha, const float* A,
           blas_int_t lda, float* B, blas_int_t ldb );

void trmm( CBLAS_ORDER order, CBLAS_SIDE side, CBLAS_UPLO uplo, CBLAS_TRANSPOSE transA,
           CBLAS_DIAG diag, blas_int_t m, blas_int_t n, double alpha, const double* A,
           blas_int_t lda, double* B, blas_int_t ldb );

void trmm( CBLAS_ORDER order, CBLAS_SIDE side, CBLAS_UPLO uplo, CBLAS_TRANSPOSE transA,
           CBLAS_DIAG diag, blas_int_t m, blas_int_t n, complex<float> alpha,
           const complex<float>* A, blas_int_t lda, complex<float>* B, blas_int_t ldb );

void trmm( CBLAS_ORDER order, CBLAS_SIDE side, CBLAS_UPLO uplo, CBLAS_TRANSPOSE transA,
           CBLAS_DIAG diag, blas_int_t m, blas_int_t n, complex<double> alpha,
           const complex<double>* A, blas_int_t lda, complex<double>* B, blas_int_t ldb );

template< typename MT1, bool SO1, typename MT2, bool SO2, typename ST >
void trmm( DenseMatrix<MT1,SO1>& B, const DenseMatrix<MT2,SO2>& A,
           CBLAS_SIDE side, CBLAS_UPLO uplo, ST alpha );

} // namespace blaze

Triangular System Solver (trsm)

The following wrapper functions provide a generic interface for the BLAS functions for solving a triangular system of equations (cblas_strsm(), cblas_dtrsm(), cblas_ctrsm(), and cblas_ztrsm()):

namespace blaze {

void trsm( CBLAS_ORDER order, CBLAS_SIDE side, CBLAS_UPLO uplo, CBLAS_TRANSPOSE transA,
           CBLAS_DIAG diag, blas_int_t m, blas_int_t n, float alpha, const float* A,
           blas_int_t lda, float* B, blas_int_t ldb );

void trsm( CBLAS_ORDER order, CBLAS_SIDE side, CBLAS_UPLO uplo, CBLAS_TRANSPOSE transA,
           CBLAS_DIAG diag, blas_int_t m, blas_int_t n, double alpha, const double* A,
           blas_int_t lda, double* B, blas_int_t ldb );

void trsm( CBLAS_ORDER order, CBLAS_SIDE side, CBLAS_UPLO uplo, CBLAS_TRANSPOSE transA,
           CBLAS_DIAG diag, blas_int_t m, blas_int_t n, complex<float> alpha,
           const complex<float>* A, blas_int_t lda, complex<float>* B, blas_int_t ldb );

void trsm( CBLAS_ORDER order, CBLAS_SIDE side, CBLAS_UPLO uplo, CBLAS_TRANSPOSE transA,
           CBLAS_DIAG diag, blas_int_t m, blas_int_t n, complex<double> alpha,
           const complex<double>* A, blas_int_t lda, complex<double>* B, blas_int_t ldb );

template< typename MT, bool SO, typename VT, bool TF, typename ST >
void trsm( const DenseMatrix<MT,SO>& A, DenseVector<VT,TF>& b,
           CBLAS_SIDE side, CBLAS_UPLO uplo, ST alpha );

template< typename MT1, bool SO1, typename MT2, bool SO2, typename ST >
void trsm( const DenseMatrix<MT1,SO1>& A, DenseMatrix<MT2,SO2>& B,
           CBLAS_SIDE side, CBLAS_UPLO uplo, ST alpha );

} // namespace blaze

Previous: Error Reporting Customization ---- Next: LAPACK Functions

Updated