-
Notifications
You must be signed in to change notification settings - Fork 3
/
Copy pathtest_linalg_eigs.f90
415 lines (308 loc) · 12.9 KB
/
test_linalg_eigs.f90
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
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
! Test eigendecomposition
module test_linalg_eig
use stdlib_linalg_interface
implicit none(type,external)
contains
!> SVD tests
subroutine test_eig(error)
logical,intent(out) :: error
real :: t0,t1
call cpu_time(t0)
call test_eig_real_s(error)
if (error) return
call test_eig_real_d(error)
if (error) return
call test_eig_real_q(error)
if (error) return
call test_eigh_real_s(error)
if (error) return
call test_eigh_real_d(error)
if (error) return
call test_eigh_real_q(error)
if (error) return
call test_eig_complex_c(error)
if (error) return
call test_eig_complex_z(error)
if (error) return
call test_eig_complex_w(error)
if (error) return
call cpu_time(t1)
print 1,1000*(t1 - t0),merge('SUCCESS','ERROR ',.not. error)
1 format('Eigenproblem tests completed in ',f9.4,' milliseconds, result=',a)
end subroutine test_eig
!> Simple real matrix eigenvalues
subroutine test_eig_real_s(error)
logical,intent(out) :: error
!> Reference solution
real(sp),parameter :: zero = 0.0_sp
real(sp),parameter :: two = 2.0_sp
real(sp),parameter :: sqrt2o2 = sqrt(two)*0.5_sp
real(sp),parameter :: tol = sqrt(epsilon(zero))
!> Local variables
type(linalg_state) :: state
real(sp) :: A(3,3),B(2,2)
complex(sp) :: lambda(3),Bvec(2,2),Bres(2,2)
!> Matrix with real eigenvalues
A = reshape([1,0,0, &
0,2,0, &
0,0,3], [3,3])
call eig(A,lambda,err=state)
error = state%error() .or. &
.not. all(aimag(lambda) == zero .and. real(lambda,kind=sp) == [1,2,3])
if (error) return
!> Matrix with complex eigenvalues
B = transpose(reshape([1,-1, &
1,1], [2,2]))
!> Expected right eigenvectors
Bres(1,1:2) = sqrt2o2
Bres(2,1) = cmplx(zero,-sqrt2o2,kind=sp)
Bres(2,2) = cmplx(zero,+sqrt2o2,kind=sp)
call eig(B,lambda,right=Bvec,err=state)
error = state%error() .or. any(abs(Bres - Bvec) > tol)
print *, bvec(1,:)
print *, bvec(2,:)
print *, bres
if (error) return
end subroutine test_eig_real_s
! Symmetric matrix eigenvalues
subroutine test_eigh_real_s(error)
logical,intent(out) :: error
!> Reference solution
real(sp),parameter :: zero = 0.0_sp
real(sp),parameter :: tol = sqrt(epsilon(zero))
real(sp),parameter :: A(4,4) = reshape([6,3,1,5, &
3,0,5,1, &
1,5,6,2, &
5,1,2,2], [4,4])
!> Local variables
real(sp) :: Amat(4,4),lambda(4),vect(4,4),Av(4,4),lv(4,4)
type(linalg_state) :: state
Amat = A
call eigh(Amat,lambda,vect,err=state)
Av = matmul(A,vect)
lv = matmul(vect,diag(lambda))
error = state%error() .or. .not. all(abs(Av - lv) < tol*abs(Av))
if (error) return
!> Test functional versions
lambda = eigvalsh(Amat)
lambda = eigvalsh(Amat,err=state)
error = state%error()
if (error) return
!> Test functional versions
Amat = A
lambda = eigvalsh(Amat,upper_a=.false.,err=state)
error = state%error()
if (error) return
end subroutine test_eigh_real_s
subroutine test_eig_real_d(error)
logical,intent(out) :: error
!> Reference solution
real(dp),parameter :: zero = 0.0_dp
real(dp),parameter :: two = 2.0_dp
real(dp),parameter :: sqrt2o2 = sqrt(two)*0.5_dp
real(dp),parameter :: tol = sqrt(epsilon(zero))
!> Local variables
type(linalg_state) :: state
real(dp) :: A(3,3),B(2,2)
complex(dp) :: lambda(3),Bvec(2,2),Bres(2,2)
!> Matrix with real eigenvalues
A = reshape([1,0,0, &
0,2,0, &
0,0,3], [3,3])
call eig(A,lambda,err=state)
error = state%error() .or. &
.not. all(aimag(lambda) == zero .and. real(lambda,kind=dp) == [1,2,3])
if (error) return
!> Matrix with complex eigenvalues
B = transpose(reshape([1,-1, &
1,1], [2,2]))
!> Expected right eigenvectors
Bres(1,1:2) = sqrt2o2
Bres(2,1) = cmplx(zero,-sqrt2o2,kind=dp)
Bres(2,2) = cmplx(zero,+sqrt2o2,kind=dp)
call eig(B,lambda,right=Bvec,err=state)
error = state%error() .or. any(abs(Bres - Bvec) > tol)
print *, bvec(1,:)
print *, bvec(2,:)
print *, bres
if (error) return
end subroutine test_eig_real_d
! Symmetric matrix eigenvalues
subroutine test_eigh_real_d(error)
logical,intent(out) :: error
!> Reference solution
real(dp),parameter :: zero = 0.0_dp
real(dp),parameter :: tol = sqrt(epsilon(zero))
real(dp),parameter :: A(4,4) = reshape([6,3,1,5, &
3,0,5,1, &
1,5,6,2, &
5,1,2,2], [4,4])
!> Local variables
real(dp) :: Amat(4,4),lambda(4),vect(4,4),Av(4,4),lv(4,4)
type(linalg_state) :: state
Amat = A
call eigh(Amat,lambda,vect,err=state)
Av = matmul(A,vect)
lv = matmul(vect,diag(lambda))
error = state%error() .or. .not. all(abs(Av - lv) < tol*abs(Av))
if (error) return
!> Test functional versions
lambda = eigvalsh(Amat)
lambda = eigvalsh(Amat,err=state)
error = state%error()
if (error) return
!> Test functional versions
Amat = A
lambda = eigvalsh(Amat,upper_a=.false.,err=state)
error = state%error()
if (error) return
end subroutine test_eigh_real_d
subroutine test_eig_real_q(error)
logical,intent(out) :: error
!> Reference solution
real(qp),parameter :: zero = 0.0_qp
real(qp),parameter :: two = 2.0_qp
real(qp),parameter :: sqrt2o2 = sqrt(two)*0.5_qp
real(qp),parameter :: tol = sqrt(epsilon(zero))
!> Local variables
type(linalg_state) :: state
real(qp) :: A(3,3),B(2,2)
complex(qp) :: lambda(3),Bvec(2,2),Bres(2,2)
!> Matrix with real eigenvalues
A = reshape([1,0,0, &
0,2,0, &
0,0,3], [3,3])
call eig(A,lambda,err=state)
error = state%error() .or. &
.not. all(aimag(lambda) == zero .and. real(lambda,kind=qp) == [1,2,3])
if (error) return
!> Matrix with complex eigenvalues
B = transpose(reshape([1,-1, &
1,1], [2,2]))
!> Expected right eigenvectors
Bres(1,1:2) = sqrt2o2
Bres(2,1) = cmplx(zero,-sqrt2o2,kind=qp)
Bres(2,2) = cmplx(zero,+sqrt2o2,kind=qp)
call eig(B,lambda,right=Bvec,err=state)
error = state%error() .or. any(abs(Bres - Bvec) > tol)
print *, bvec(1,:)
print *, bvec(2,:)
print *, bres
if (error) return
end subroutine test_eig_real_q
! Symmetric matrix eigenvalues
subroutine test_eigh_real_q(error)
logical,intent(out) :: error
!> Reference solution
real(qp),parameter :: zero = 0.0_qp
real(qp),parameter :: tol = sqrt(epsilon(zero))
real(qp),parameter :: A(4,4) = reshape([6,3,1,5, &
3,0,5,1, &
1,5,6,2, &
5,1,2,2], [4,4])
!> Local variables
real(qp) :: Amat(4,4),lambda(4),vect(4,4),Av(4,4),lv(4,4)
type(linalg_state) :: state
Amat = A
call eigh(Amat,lambda,vect,err=state)
Av = matmul(A,vect)
lv = matmul(vect,diag(lambda))
error = state%error() .or. .not. all(abs(Av - lv) < tol*abs(Av))
if (error) return
!> Test functional versions
lambda = eigvalsh(Amat)
lambda = eigvalsh(Amat,err=state)
error = state%error()
if (error) return
!> Test functional versions
Amat = A
lambda = eigvalsh(Amat,upper_a=.false.,err=state)
error = state%error()
if (error) return
end subroutine test_eigh_real_q
!> Simple complex matrix eigenvalues
subroutine test_eig_complex_c(error)
logical,intent(out) :: error
!> Reference solution
real(sp),parameter :: zero = 0.0_sp
real(sp),parameter :: two = 2.0_sp
real(sp),parameter :: sqrt2o2 = sqrt(two)*0.5_sp
real(sp),parameter :: tol = sqrt(epsilon(zero))
complex(sp),parameter :: cone = (1.0_sp,0.0_sp)
complex(sp),parameter :: cimg = (0.0_sp,1.0_sp)
complex(sp),parameter :: czero = (0.0_sp,0.0_sp)
!> Local vaciables
type(linalg_state) :: state
complex(sp) :: A(2,2),lambda(2),Avec(2,2),Ares(2,2),lres(2)
!> Matcix with real eigenvalues
A = transpose(reshape([cone,cimg, &
-cimg,cone], [2,2]))
call eig(A,lambda,right=Avec,err=state)
!> Expected eigenvalues and eigenvectors
lres(1) = two
lres(2) = zero
!> Eigenvectors may vary: do not use for error
Ares(1,1) = cmplx(zero,sqrt2o2,kind=sp)
Ares(1,2) = cmplx(sqrt2o2,zero,kind=sp)
Ares(2,1) = cmplx(sqrt2o2,zero,kind=sp)
Ares(2,2) = cmplx(zero,sqrt2o2,kind=sp)
error = state%error() .or. any(abs(lambda - lres) > tol)
if (error) return
end subroutine test_eig_complex_c
subroutine test_eig_complex_z(error)
logical,intent(out) :: error
!> Reference solution
real(dp),parameter :: zero = 0.0_dp
real(dp),parameter :: two = 2.0_dp
real(dp),parameter :: sqrt2o2 = sqrt(two)*0.5_dp
real(dp),parameter :: tol = sqrt(epsilon(zero))
complex(dp),parameter :: cone = (1.0_dp,0.0_dp)
complex(dp),parameter :: cimg = (0.0_dp,1.0_dp)
complex(dp),parameter :: czero = (0.0_dp,0.0_dp)
!> Local vaciables
type(linalg_state) :: state
complex(dp) :: A(2,2),lambda(2),Avec(2,2),Ares(2,2),lres(2)
!> Matcix with real eigenvalues
A = transpose(reshape([cone,cimg, &
-cimg,cone], [2,2]))
call eig(A,lambda,right=Avec,err=state)
!> Expected eigenvalues and eigenvectors
lres(1) = two
lres(2) = zero
!> Eigenvectors may vary: do not use for error
Ares(1,1) = cmplx(zero,sqrt2o2,kind=dp)
Ares(1,2) = cmplx(sqrt2o2,zero,kind=dp)
Ares(2,1) = cmplx(sqrt2o2,zero,kind=dp)
Ares(2,2) = cmplx(zero,sqrt2o2,kind=dp)
error = state%error() .or. any(abs(lambda - lres) > tol)
if (error) return
end subroutine test_eig_complex_z
subroutine test_eig_complex_w(error)
logical,intent(out) :: error
!> Reference solution
real(qp),parameter :: zero = 0.0_qp
real(qp),parameter :: two = 2.0_qp
real(qp),parameter :: sqrt2o2 = sqrt(two)*0.5_qp
real(qp),parameter :: tol = sqrt(epsilon(zero))
complex(qp),parameter :: cone = (1.0_qp,0.0_qp)
complex(qp),parameter :: cimg = (0.0_qp,1.0_qp)
complex(qp),parameter :: czero = (0.0_qp,0.0_qp)
!> Local vaciables
type(linalg_state) :: state
complex(qp) :: A(2,2),lambda(2),Avec(2,2),Ares(2,2),lres(2)
!> Matcix with real eigenvalues
A = transpose(reshape([cone,cimg, &
-cimg,cone], [2,2]))
call eig(A,lambda,right=Avec,err=state)
!> Expected eigenvalues and eigenvectors
lres(1) = two
lres(2) = zero
!> Eigenvectors may vary: do not use for error
Ares(1,1) = cmplx(zero,sqrt2o2,kind=qp)
Ares(1,2) = cmplx(sqrt2o2,zero,kind=qp)
Ares(2,1) = cmplx(sqrt2o2,zero,kind=qp)
Ares(2,2) = cmplx(zero,sqrt2o2,kind=qp)
error = state%error() .or. any(abs(lambda - lres) > tol)
if (error) return
end subroutine test_eig_complex_w
end module test_linalg_eig