mirror of
https://github.com/TREX-CoE/qmckl.git
synced 2025-01-03 18:16:28 +01:00
Tests pass for qmckl_dgemm. #32.
This commit is contained in:
parent
cf3550b6b7
commit
eaede28a73
@ -90,13 +90,13 @@ integer function qmckl_dgemm_f(context, TransA, TransB, m, n, k, alpha, A, LDA,
|
||||
integer*8 , intent(in) :: m, n, k
|
||||
real*8 , intent(in) :: alpha, beta
|
||||
integer*8 , intent(in) :: lda
|
||||
real*8 , intent(in) :: A(m,n)
|
||||
real*8 , intent(in) :: A(m,k)
|
||||
integer*8 , intent(in) :: ldb
|
||||
real*8 , intent(in) :: B(n,k)
|
||||
real*8 , intent(in) :: B(k,n)
|
||||
integer*8 , intent(in) :: ldc
|
||||
real*8 , intent(out) :: C(m,n)
|
||||
|
||||
integer*8 :: i,j
|
||||
integer*8 :: i,j,l
|
||||
|
||||
info = QMCKL_SUCCESS
|
||||
|
||||
@ -120,17 +120,17 @@ integer function qmckl_dgemm_f(context, TransA, TransB, m, n, k, alpha, A, LDA,
|
||||
return
|
||||
endif
|
||||
|
||||
if (LDA < m) then
|
||||
if (LDA .ne. m) then
|
||||
info = QMCKL_INVALID_ARG_9
|
||||
return
|
||||
endif
|
||||
|
||||
if (LDB < n) then
|
||||
if (LDB .ne. k) then
|
||||
info = QMCKL_INVALID_ARG_10
|
||||
return
|
||||
endif
|
||||
|
||||
if (LDB < n) then
|
||||
if (LDC .ne. m) then
|
||||
info = QMCKL_INVALID_ARG_13
|
||||
return
|
||||
endif
|
||||
@ -153,8 +153,8 @@ end function qmckl_dgemm_f
|
||||
implicit none
|
||||
|
||||
integer (c_int64_t) , intent(in) , value :: context
|
||||
logical*8 , intent(in) , value :: TransA
|
||||
logical*8 , intent(in) , value :: TransB
|
||||
logical , intent(in) , value :: TransA
|
||||
logical , intent(in) , value :: TransB
|
||||
integer (c_int64_t) , intent(in) , value :: m
|
||||
integer (c_int64_t) , intent(in) , value :: n
|
||||
integer (c_int64_t) , intent(in) , value :: k
|
||||
@ -188,17 +188,17 @@ end function qmckl_dgemm_f
|
||||
implicit none
|
||||
|
||||
integer (c_int64_t) , intent(in) , value :: context
|
||||
logical*8 , intent(in) , value :: TransA
|
||||
logical*8 , intent(in) , value :: TransB
|
||||
logical , intent(in) , value :: TransA
|
||||
logical , intent(in) , value :: TransB
|
||||
integer (c_int64_t) , intent(in) , value :: m
|
||||
integer (c_int64_t) , intent(in) , value :: n
|
||||
integer (c_int64_t) , intent(in) , value :: k
|
||||
real (c_double ) , intent(in) :: alpha
|
||||
real (c_double ) , intent(in) , value :: alpha
|
||||
integer (c_int64_t) , intent(in) , value :: lda
|
||||
real (c_double ) , intent(in) :: A(lda,*)
|
||||
integer (c_int64_t) , intent(in) , value :: ldb
|
||||
real (c_double ) , intent(in) :: B(ldb,*)
|
||||
real (c_double ) , intent(in) :: beta
|
||||
real (c_double ) , intent(in) , value :: beta
|
||||
integer (c_int64_t) , intent(in) , value :: ldc
|
||||
real (c_double ) , intent(out) :: C(ldc,*)
|
||||
|
||||
@ -216,7 +216,7 @@ integer(qmckl_exit_code) function test_qmckl_dgemm(context) bind(C)
|
||||
double precision, allocatable :: A(:,:), B(:,:), C(:,:), D(:,:)
|
||||
integer*8 :: m, n, k, LDA, LDB, LDC
|
||||
integer*8 :: i,j,l
|
||||
logical*8 :: TransA, TransB
|
||||
logical :: TransA, TransB
|
||||
double precision :: x, alpha, beta
|
||||
|
||||
TransA = .False.
|
||||
@ -234,14 +234,16 @@ integer(qmckl_exit_code) function test_qmckl_dgemm(context) bind(C)
|
||||
B = 0.d0
|
||||
C = 0.d0
|
||||
D = 0.d0
|
||||
do j=1,m
|
||||
do i=1,k
|
||||
alpha = 1.0d0
|
||||
beta = 0.0d0
|
||||
do j=1,k
|
||||
do i=1,m
|
||||
A(i,j) = -10.d0 + dble(i+j)
|
||||
end do
|
||||
end do
|
||||
|
||||
do j=1,k
|
||||
do i=1,n
|
||||
do j=1,n
|
||||
do i=1,k
|
||||
B(i,j) = -10.d0 + dble(i+j)
|
||||
end do
|
||||
end do
|
||||
@ -253,21 +255,20 @@ integer(qmckl_exit_code) function test_qmckl_dgemm(context) bind(C)
|
||||
test_qmckl_dgemm = QMCKL_FAILURE
|
||||
|
||||
x = 0.d0
|
||||
do j=1,m
|
||||
do i=l,n
|
||||
do j=1,n
|
||||
do i=1,m
|
||||
do l=1,k
|
||||
D(i,j) = D(i,j) + A(i,k)*B(k,j)
|
||||
D(i,j) = D(i,j) + A(i,l)*B(l,j)
|
||||
end do
|
||||
x = x + (D(i,j) - C(i,j))**2
|
||||
end do
|
||||
end do
|
||||
print *,"DABS(X)= ",dabs(x)
|
||||
|
||||
if (dabs(x) <= 1.d-15) then
|
||||
test_qmckl_dgemm = QMCKL_SUCCESS
|
||||
endif
|
||||
|
||||
deallocate(A,B)
|
||||
deallocate(A,B,C,D)
|
||||
end function test_qmckl_dgemm
|
||||
#+end_src
|
||||
|
||||
|
Loading…
Reference in New Issue
Block a user