/
dsbevd_example.f90
executable file
·83 lines (69 loc) · 2.39 KB
/
dsbevd_example.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
Program dsbevd_example
! DSBEVD Example Program Text
! Copyright (c) 2018, Numerical Algorithms Group (NAG Ltd.)
! For licence see
! https://github.com/numericalalgorithmsgroup/LAPACK_Examples/blob/master/LICENCE.md
! .. Use Statements ..
Use lapack_example_aux, Only: nagf_file_print_matrix_real_gen
Use lapack_interfaces, Only: dsbevd
Use lapack_precision, Only: dp
! .. Implicit None Statement ..
Implicit None
! .. Parameters ..
Integer, Parameter :: nin = 5, nout = 6
! .. Local Scalars ..
Integer :: i, ifail, info, j, kd, ldab, ldz, liwork, lwork, n
Character (1) :: job, uplo
! .. Local Arrays ..
Real (Kind=dp), Allocatable :: ab(:, :), w(:), work(:), z(:, :)
Integer, Allocatable :: iwork(:)
! .. Intrinsic Procedures ..
Intrinsic :: max, min
! .. Executable Statements ..
Write (nout, *) 'DSBEVD Example Program Results'
! Skip heading in data file
Read (nin, *)
Read (nin, *) n, kd
ldab = kd + 1
ldz = n
liwork = 5*n + 3
lwork = 2*n*n + 5*n + 1
Allocate (ab(ldab,n), w(n), work(lwork), z(ldz,n), iwork(liwork))
! Read A from data file
Read (nin, *) uplo
If (uplo=='U') Then
Do i = 1, n
Read (nin, *)(ab(kd+1+i-j,j), j=i, min(n,i+kd))
End Do
Else If (uplo=='L') Then
Do i = 1, n
Read (nin, *)(ab(1+i-j,j), j=max(1,i-kd), i)
End Do
End If
Read (nin, *) job
! Calculate all the eigenvalues and eigenvectors of A
Call dsbevd(job, uplo, n, kd, ab, ldab, w, z, ldz, work, lwork, iwork, &
liwork, info)
Write (nout, *)
If (info>0) Then
Write (nout, *) 'Failure to converge.'
Else
! Print eigenvalues and eigenvectors
Write (nout, *) 'Eigenvalues'
Write (nout, 100) w(1:n)
Write (nout, *)
Flush (nout)
! Standardize the eigenvectors so that first elements are non-negative.
Do i = 1, n
If (z(1,i)<0.0_dp) Then
z(1:n, i) = -z(1:n, i)
End If
End Do
! ifail: behaviour on error exit
! =0 for hard exit, =1 for quiet-soft, =-1 for noisy-soft
ifail = 0
Call nagf_file_print_matrix_real_gen('General', ' ', n, n, z, ldz, &
'Eigenvectors', ifail)
End If
100 Format (3X, (8F8.4))
End Program