-
Notifications
You must be signed in to change notification settings - Fork 185
/
Copy pathtest_blas_lapack.fypp
218 lines (166 loc) · 6.57 KB
/
test_blas_lapack.fypp
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
#:include "common.fypp"
#:set RC_KINDS_TYPES = REAL_KINDS_TYPES + CMPLX_KINDS_TYPES
module test_blas_lapack
use testdrive, only : new_unittest, unittest_type, error_type, check, skip_test
use stdlib_kinds, only: sp, dp, xdp, qp, int8, int16, int32, int64
use stdlib_linalg, only: eye
use stdlib_linalg_blas
use stdlib_linalg_lapack
implicit none
real(sp), parameter :: sptol = 1000 * epsilon(1._sp)
real(dp), parameter :: dptol = 1000 * epsilon(1._dp)
#:if WITH_QP
real(qp), parameter :: qptol = 1000 * epsilon(1._qp)
#:endif
contains
!> Collect all exported unit tests
subroutine collect_blas_lapack(testsuite)
!> Collection of tests
type(unittest_type), allocatable, intent(out) :: testsuite(:)
testsuite = [ &
#:for k1, t1 in REAL_KINDS_TYPES
new_unittest("test_gemv${t1[0]}$${k1}$", test_gemv${t1[0]}$${k1}$), &
new_unittest("test_getri${t1[0]}$${k1}$", test_getri${t1[0]}$${k1}$), &
#:endfor
new_unittest("test_idamax", test_idamax), &
new_unittest("test_external_blas",external_blas_test), &
new_unittest("test_external_lapack",external_lapack_test) &
]
end subroutine collect_blas_lapack
#:for k1, t1 in REAL_KINDS_TYPES
subroutine test_gemv${t1[0]}$${k1}$(error)
!> Error handling
type(error_type), allocatable, intent(out) :: error
#:if k1=="xdp"
call skip_test(error, "Extended precision is not enabled")
#:else
${t1}$ :: A(3,3),x(3),y(3),ylap(3),yintr(3),alpha,beta
call random_number(alpha)
call random_number(beta)
call random_number(A)
call random_number(x)
call random_number(y)
ylap = y
call gemv('No transpose',size(A,1),size(A,2),alpha,A,size(A,1),x,1,beta,ylap,1)
yintr = alpha*matmul(A,x)+beta*y
call check(error, sum(abs(ylap - yintr)) < sptol, &
"blas vs. intrinsics axpy: sum() < sptol failed")
if (allocated(error)) return
#:endif
end subroutine test_gemv${t1[0]}$${k1}$
! Find matrix inverse from LU decomposition
subroutine test_getri${t1[0]}$${k1}$(error)
!> Error handling
type(error_type), allocatable, intent(out) :: error
#:if k1=="xdp"
call skip_test(error, "Extended precision is not enabled")
#:else
integer(ilp), parameter :: n = 3
${t1}$ :: A(n,n)
${t1}$,allocatable :: work(:)
integer(ilp) :: ipiv(n),info,lwork,nb
A = eye(n)
! Factorize matrix (overwrite result)
call getrf(size(A,1),size(A,2),A,size(A,1),ipiv,info)
call check(error, info==0, "lapack getrf returned info/=0")
if (allocated(error)) return
! Get optimal worksize (returned in work(1)) (apply 2% safety parameter)
nb = stdlib_ilaenv(1,'${t1[0]}$getri',' ',n,-1,-1,-1)
lwork = nint(1.02*n*nb,kind=ilp)
allocate (work(lwork))
! Invert matrix
call getri(n,a,n,ipiv,work,lwork,info)
call check(error, info==0, "lapack getri returned info/=0")
if (allocated(error)) return
call check(error, sum(abs(A - eye(3))) < sptol, &
"lapack eye inversion: tolerance check failed")
if (allocated(error)) return
#:endif
end subroutine test_getri${t1[0]}$${k1}$
#:endfor
! Return
subroutine test_idamax(error)
!> Error handling
type(error_type), allocatable, intent(out) :: error
integer(ilp), parameter :: n = 5
integer(ilp) :: imax
real(dp) :: x(n)
x = [1,2,3,4,5]
imax = stdlib_idamax(n,x,1)
call check(error, imax==5, "blas idamax returned wrong location")
end subroutine test_idamax
!> Test availability of the external BLAS interface
subroutine external_blas_test(error)
!> Error handling
type(error_type), allocatable, intent(out) :: error
#ifdef STDLIB_EXTERNAL_BLAS
interface
subroutine saxpy(n,sa,sx,incx,sy,incy)
import sp,ilp
implicit none(type,external)
real(sp), intent(in) :: sa,sx(*)
integer(ilp), intent(in) :: incx,incy,n
real(sp), intent(inout) :: sy(*)
end subroutine saxpy
end interface
integer(ilp), parameter :: n = 5, inc=1
real(sp) :: a,x(n),y(n)
x = 1.0_sp
y = 2.0_sp
a = 3.0_sp
call saxpy(n,a,x,inc,y,inc)
call check(error, all(abs(y-5.0_sp)<sqrt(epsilon(0.0_sp))), "saxpy: check result")
if (allocated(error)) return
#else
call skip_test(error, "Not using an external BLAS")
#endif
end subroutine external_blas_test
!> Test availability of the external BLAS interface
subroutine external_lapack_test(error)
!> Error handling
type(error_type), allocatable, intent(out) :: error
#ifdef STDLIB_EXTERNAL_LAPACK
interface
subroutine dgetrf( m, n, a, lda, ipiv, info )
import dp,ilp
implicit none(type,external)
integer(ilp), intent(out) :: info,ipiv(*)
integer(ilp), intent(in) :: lda,m,n
real(dp), intent(inout) :: a(lda,*)
end subroutine dgetrf
end interface
integer(ilp), parameter :: n = 3
real(dp) :: A(n,n)
integer(ilp) :: ipiv(n),info
A = eye(n)
info = 123
! Factorize matrix
call dgetrf(n,n,A,n,ipiv,info)
call check(error, info==0, "dgetrf: check result")
if (allocated(error)) return
#else
call skip_test(error, "Not using an external LAPACK")
#endif
end subroutine external_lapack_test
end module test_blas_lapack
program tester
use, intrinsic :: iso_fortran_env, only : error_unit
use testdrive, only : run_testsuite, new_testsuite, testsuite_type
use test_blas_lapack, only : collect_blas_lapack
implicit none
integer :: stat, is
type(testsuite_type), allocatable :: testsuites(:)
character(len=*), parameter :: fmt = '("#", *(1x, a))'
stat = 0
testsuites = [ &
new_testsuite("blas_lapack", collect_blas_lapack) &
]
do is = 1, size(testsuites)
write(error_unit, fmt) "Testing:", testsuites(is)%name
call run_testsuite(testsuites(is)%collect, error_unit, stat)
end do
if (stat > 0) then
write(error_unit, '(i0, 1x, a)') stat, "test(s) failed!"
error stop
end if
end program