BLAS functions
Table of Contents
1 Matrix operations
1.1 qmckl_dgemm
Matrix multiply l\(C_{ij} = \beta C_{ij} + \alpha \sum_{k} A_{ik} \cdot B_{kj}\) using Fortran matmul
function.
TODO: Add description about the external library dependence.
qmcklcontext | context | in | Global state |
bool | TransA | in | Number of rows of the input matrix |
bool | TransB | in | Number of rows of the input matrix |
int64t | m | in | Number of rows of the input matrix |
int64t | n | in | Number of columns of the input matrix |
int64t | k | in | Number of columns of the input matrix |
double | alpha | in | Number of columns of the input matrix |
double | A[][lda] | in | Array containing the \(m \times n\) matrix \(A\) |
int64t | lda | in | Leading dimension of array A |
double | B[][ldb] | in | Array containing the \(n \times m\) matrix \(B\) |
int64t | ldb | in | Leading dimension of array B |
double | beta | in | Array containing the \(n \times m\) matrix \(B\) |
double | C[][ldc] | out | Array containing the \(n \times m\) matrix \(B\) |
int64t | ldc | in | Leading dimension of array B |
1.1.1 Requirements
context
is notQMCKL_NULL_CONTEXT
m > 0
n > 0
k > 0
lda >= m
ldb >= n
ldc >= n
A
is allocated with at least \(m \times k \times 8\) bytesB
is allocated with at least \(k \times n \times 8\) bytesC
is allocated with at least \(m \times n \times 8\) bytes
1.1.2 C header
qmckl_exit_code qmckl_dgemm ( const qmckl_context context, const bool TransA, const bool TransB, const int64_t m, const int64_t n, const int64_t k, const double alpha, const double* A, const int64_t lda, const double* B, const int64_t ldb, const double beta, double* const C, const int64_t ldc );
1.1.3 Source
integer function qmckl_dgemm_f(context, TransA, TransB, m, n, k, alpha, A, LDA, B, LDB, beta, C, LDC) & result(info) use qmckl implicit none integer(qmckl_context) , intent(in) :: context logical*8 , intent(in) :: TransA, TransB integer*8 , intent(in) :: m, n, k real*8 , intent(in) :: alpha, beta integer*8 , intent(in) :: lda real*8 , intent(in) :: A(m,k) integer*8 , intent(in) :: ldb real*8 , intent(in) :: B(k,n) integer*8 , intent(in) :: ldc real*8 , intent(out) :: C(m,n) real*8, allocatable :: AT(:,:), BT(:,:), CT(:,:) integer*8 :: i,j,l, LDA_2, LDB_2 info = QMCKL_SUCCESS if (TransA) then allocate(AT(k,m)) do i = 1, m do j = 1, k AT(j,i) = A(i,j) end do end do LDA_2 = M else LDA_2 = LDA endif if (TransB) then allocate(BT(n,k)) do i = 1, k do j = 1, n BT(j,i) = B(i,j) end do end do LDB_2 = K else LDB_2 = LDB endif if (context == QMCKL_NULL_CONTEXT) then info = QMCKL_INVALID_CONTEXT return endif if (m <= 0_8) then info = QMCKL_INVALID_ARG_4 return endif if (n <= 0_8) then info = QMCKL_INVALID_ARG_5 return endif if (k <= 0_8) then info = QMCKL_INVALID_ARG_6 return endif if (LDA_2 .ne. m) then info = QMCKL_INVALID_ARG_9 return endif if (LDB_2 .ne. k) then info = QMCKL_INVALID_ARG_10 return endif if (LDC .ne. m) then info = QMCKL_INVALID_ARG_13 return endif if (TransA) then C = matmul(AT,B) else if (TransB) then C = matmul(A,BT) else C = matmul(A,B) endif end function qmckl_dgemm_f