-
Notifications
You must be signed in to change notification settings - Fork 2
/
Copy pathdsyevj3.f
167 lines (153 loc) · 5.65 KB
/
dsyevj3.f
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
!DEC$ IF .NOT. DEFINED (dsyevj3_F)
!DEC$ DEFINE dsyevj3_F
* ----------------------------------------------------------------------------
* Numerical diagonalization of 3x3 matrcies
* Copyright (C) 2006 Joachim Kopp
* [https://www.mpi-hd.mpg.de/personalhomes/globes/3x3/index.html]
* ----------------------------------------------------------------------------
* This library is free software; you can redistribute it and/or
* modify it under the terms of the GNU Lesser General Public
* License as published by the Free Software Foundation; either
* version 2.1 of the License, or (at your option) any later version.
*
* This library is distributed in the hope that it will be useful,
* but WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
* Lesser General Public License for more details.
*
* You should have received a copy of the GNU Lesser General Public
* License along with this library; if not, write to the Free Software
* Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
* ----------------------------------------------------------------------------
* ----------------------------------------------------------------------------
SUBROUTINE dsyevj3(A_in, Q, W, failedIt)
* ----------------------------------------------------------------------------
* Calculates the eigenvalues and normalized eigenvectors of a symmetric 3x3
* matrix A using the Jacobi algorithm.
* The upper triangular part of A is destroyed during the calculation,
* the diagonal elements are read but not destroyed, and the lower
* triangular elements are not referenced at all.
* ----------------------------------------------------------------------------
* Parameters:
* A: The symmetric input matrix
* Q: Storage buffer for eigenvectors
* W: Storage buffer for eigenvalues
* ----------------------------------------------------------------------------
* .. Arguments ..
DOUBLE PRECISION, intent(in) :: A_in(3,3)
DOUBLE PRECISION, intent(out) :: Q(3,3)
DOUBLE PRECISION, intent(out) :: W(3)
logical, intent(out) :: failedIt
* .. Parameters ..
INTEGER N
PARAMETER ( N = 3 )
* .. Local Variables ..
DOUBLE PRECISION A(3,3)
DOUBLE PRECISION SD, SO
DOUBLE PRECISION S, C, T
DOUBLE PRECISION G, H, Z, THETA
DOUBLE PRECISION THRESH
INTEGER I, X, Y, R
c
failedIt = .false.
do X=1,N
do Y=1,N
A(X,Y) = A_in(X,Y)
enddo
enddo
c
* Initialize Q to the identitity matrix
* --- This loop can be omitted if only the eigenvalues are desired ---
DO 10 X = 1, N
Q(X,X) = 1.0D0
DO 11, Y = 1, X-1
Q(X, Y) = 0.0D0
Q(Y, X) = 0.0D0
11 CONTINUE
10 CONTINUE
* Initialize W to diag(A)
DO 20 X = 1, N
W(X) = A(X, X)
20 CONTINUE
* Calculate SQR(tr(A))
SD = 0.0D0
DO 30 X = 1, N
SD = SD + ABS(W(X))
30 CONTINUE
SD = SD**2
* Main iteration loop
DO 40 I = 1, 50
* Test for convergence
SO = 0.0D0
DO 50 X = 1, N
DO 51 Y = X+1, N
SO = SO + ABS(A(X, Y))
51 CONTINUE
50 CONTINUE
IF (SO .EQ. 0.0D0) THEN
RETURN
END IF
IF (I .LT. 4) THEN
THRESH = 0.2D0 * SO / N**2
ELSE
THRESH = 0.0D0
END IF
* Do sweep
DO 60 X = 1, N
DO 61 Y = X+1, N
G = 100.0D0 * ( ABS(A(X, Y)) )
IF ( I .GT. 4 .AND. ABS(W(X)) + G .EQ. ABS(W(X))
$ .AND. ABS(W(Y)) + G .EQ. ABS(W(Y)) ) THEN
A(X, Y) = 0.0D0
ELSE IF (ABS(A(X, Y)) .GT. THRESH) THEN
* Calculate Jacobi transformation
H = W(Y) - W(X)
IF ( ABS(H) + G .EQ. ABS(H) ) THEN
T = A(X, Y) / H
ELSE
THETA = 0.5D0 * H / A(X, Y)
IF (THETA .LT. 0.0D0) THEN
T = -1.0D0 / (SQRT(1.0D0 + THETA**2) - THETA)
ELSE
T = 1.0D0 / (SQRT(1.0D0 + THETA**2) + THETA)
END IF
END IF
C = 1.0D0 / SQRT( 1.0D0 + T**2 )
S = T * C
Z = T * A(X, Y)
* Apply Jacobi transformation
A(X, Y) = 0.0D0
W(X) = W(X) - Z
W(Y) = W(Y) + Z
DO 70 R = 1, X-1
T = A(R, X)
A(R, X) = C * T - S * A(R, Y)
A(R, Y) = S * T + C * A(R, Y)
70 CONTINUE
DO 80, R = X+1, Y-1
T = A(X, R)
A(X, R) = C * T - S * A(R, Y)
A(R, Y) = S * T + C * A(R, Y)
80 CONTINUE
DO 90, R = Y+1, N
T = A(X, R)
A(X, R) = C * T - S * A(Y, R)
A(Y, R) = S * T + C * A(Y, R)
90 CONTINUE
* Update eigenvectors
* --- This loop can be omitted if only the eigenvalues are desired ---
DO 100, R = 1, N
T = Q(R, X)
Q(R, X) = C * T - S * Q(R, Y)
Q(R, Y) = S * T + C * Q(R, Y)
100 CONTINUE
END IF
61 CONTINUE
60 CONTINUE
40 CONTINUE
! call usermsg('DSYEVJ3<< No convergence for EV/EW.')
! write(*,*) "DSYEVJ3<< Input A_in=",A_in
failedIt = .true.
END SUBROUTINE
* End of subroutine dsyevj3
!DEC$ ENDIF