root / src / Mat_util.f90
Historique | Voir | Annoter | Télécharger (8,37 ko)
1 |
!---------------------------------------------------------------------- |
---|---|
2 |
! Copyright 2003-2014 Ecole Normale Supérieure de Lyon, |
3 |
! Centre National de la Recherche Scientifique, |
4 |
! Université Claude Bernard Lyon 1. All rights reserved. |
5 |
! |
6 |
! This work is registered with the Agency for the Protection of Programs |
7 |
! as IDDN.FR.001.100009.000.S.P.2014.000.30625 |
8 |
! |
9 |
! Authors: P. Fleurat-Lessard, P. Dayal |
10 |
! Contact: optnpath@gmail.com |
11 |
! |
12 |
! This file is part of "Opt'n Path". |
13 |
! |
14 |
! "Opt'n Path" is free software: you can redistribute it and/or modify |
15 |
! it under the terms of the GNU Affero General Public License as |
16 |
! published by the Free Software Foundation, either version 3 of the License, |
17 |
! or (at your option) any later version. |
18 |
! |
19 |
! "Opt'n Path" is distributed in the hope that it will be useful, |
20 |
! but WITHOUT ANY WARRANTY; without even the implied warranty of |
21 |
! |
22 |
! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the |
23 |
! GNU Affero General Public License for more details. |
24 |
! |
25 |
! You should have received a copy of the GNU Affero General Public License |
26 |
! along with "Opt'n Path". If not, see <http://www.gnu.org/licenses/>. |
27 |
! |
28 |
! Contact The Office of Technology Licensing, valorisation@ens-lyon.fr, |
29 |
! for commercial licensing opportunities. |
30 |
!---------------------------------------------------------------------- |
31 |
|
32 |
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
33 |
! |
34 |
SUBROUTINE GenInv(N,A,InvA,NReal) |
35 |
!!!!!!!!!!!!!!!! |
36 |
! |
37 |
! This subroutines calculates the generalized inverse of a matrix |
38 |
! It first diagonalize the matrix A, then inverse all non-zero |
39 |
! eigenvalues, and forms the InvA matrix using these new eigenvalues |
40 |
! |
41 |
! Input: |
42 |
! N : dimension of A |
43 |
! NReal :Actual dimension of A |
44 |
! A(N,N) : Initial Matrix, stored in A(Nreal,Nreal) |
45 |
! |
46 |
! Output: |
47 |
! InvA(N,N) : Inversed Matrix, stored in a (Nreal,NReal) matrix |
48 |
! |
49 |
!!!!!!!!!!!!!!!!!!!!!!! |
50 |
|
51 |
Use Vartypes |
52 |
IMPLICIT NONE |
53 |
|
54 |
INTEGER(KINT), INTENT(IN) :: N,Nreal |
55 |
REAL(KREAL), INTENT(IN) :: A(NReal,NReal) |
56 |
REAL(KREAL), INTENT(OUT) :: InvA(NReal,NReal) |
57 |
! |
58 |
|
59 |
INTEGER(KINT) :: I,J,K |
60 |
REAL(KREAL), ALLOCATABLE :: EigVec(:,:) ! (Nreal,Nreal) |
61 |
REAL(KREAL), ALLOCATABLE :: EigVal(:) ! (Nreal) |
62 |
REAL(KREAL), ALLOCATABLE :: ATmp(:,:) ! (NReal,Nreal) |
63 |
REAL(KREAL) :: ss |
64 |
! |
65 |
REAL(KREAL), PARAMETER :: eps=1e-12 |
66 |
|
67 |
ALLOCATE(EigVec(Nreal,Nreal), EigVal(Nreal),ATmp(NReal,NReal)) |
68 |
! A will be destroyed in Jacobi so we save it |
69 |
ATmp=A |
70 |
CALL JAcobi(ATmp,N,EigVal,EigVec,NReal) |
71 |
|
72 |
DO I=1,N |
73 |
IF (abs(EigVal(I)).GT.eps) EigVal(I)=1.d0/EigVal(I) |
74 |
END DO |
75 |
|
76 |
InvA=0.d0 |
77 |
do k = 1, n |
78 |
do j = 1, n |
79 |
ss = eigval(k) * eigvec(j,k) |
80 |
do i = 1, j |
81 |
InvA(i,j) = InvA(i,j) + ss * eigvec(i,k) |
82 |
end do |
83 |
end do |
84 |
end do |
85 |
do j = 1, n |
86 |
do i = 1, j-1 |
87 |
InvA(j,i) = InvA(i,j) |
88 |
end do |
89 |
end do |
90 |
|
91 |
|
92 |
DEALLOCATE(EigVec, EigVal,ATmp) |
93 |
END SUBROUTINE GenInv |
94 |
|
95 |
!============================================================ |
96 |
! |
97 |
! ++ Matrix diagonalization Using jacobi |
98 |
! Works only for symmetric matrices |
99 |
! EigvenVectors : V(i,i) |
100 |
! Eigenvalues : D(i) |
101 |
! PFL 30/05/03 |
102 |
! This version uses packed matrices. |
103 |
! it unpacks them before calling Jacobi ! |
104 |
! we have AP(i + (j-1)*j/2) = A(i,j) for 1<=i<=j; |
105 |
! |
106 |
!============================================================ |
107 |
! |
108 |
SUBROUTINE JacPacked(N,AP,D,V,nreal) |
109 |
|
110 |
Use Vartypes |
111 |
IMPLICIT NONE |
112 |
|
113 |
INTEGER(KINT), INTENT(IN) :: N,NREAL |
114 |
REAL(KREAL) :: AP(N*(N+1)/2) |
115 |
REAL(KREAL), ALLOCATABLE :: A(:,:) |
116 |
REAL(KREAL) :: V(Nreal,Nreal),D(Nreal) |
117 |
INTEGER(KINT) :: i,j,nn |
118 |
|
119 |
allocate(A(nreal,nreal)) |
120 |
nn=n*(n+1)/2 |
121 |
! WRITE(*,*) 'Jacpa 0' |
122 |
! WRITE(*,'(I3,10(1X,F15.6))') n,(AP(i),i=1,min(nn,5)) |
123 |
! WRITE(*,*) 'Jacpa 0' |
124 |
do j=1,n |
125 |
do i=1,j |
126 |
! WRITE(*,*) i,j |
127 |
A(i,J)=AP(i + (j-1)*j/2) |
128 |
A(j,I)=A(i,J) |
129 |
end do |
130 |
end do |
131 |
! do j=1,n |
132 |
! WRITE(*,'(10(1X,F15.6))') (A(i,J),i=1,min(n,5)) |
133 |
! end do |
134 |
|
135 |
! WRITE(*,*) 'Jacpa 1' |
136 |
call Jacobi(A,n,D,V,Nreal) |
137 |
! WRITE(*,*) 'Jacpa 2' |
138 |
! DO i=1,n |
139 |
! WRITE(*,'(1X,I3,10(1X,F15.6))') i,D(i),(V(j,i),j=1,min(n,5)) |
140 |
! end do |
141 |
deallocate(a) |
142 |
|
143 |
end SUBROUTINE JacPacked |
144 |
|
145 |
|
146 |
|
147 |
! |
148 |
!============================================================ |
149 |
! |
150 |
! ++ Matrix diagonalization Using jacobi |
151 |
! Works only for symmetric matrices |
152 |
! EigvenVectors : V |
153 |
! Eigenvalues : D |
154 |
! |
155 |
!============================================================ |
156 |
! |
157 |
SUBROUTINE JACOBI(A,N,D,V,Nreal) |
158 |
|
159 |
!!!!!!!!!!!!!!!! |
160 |
! |
161 |
! Input: |
162 |
! N : Dimension of A |
163 |
! NReal : Actual dimensions of A, D and V. |
164 |
! |
165 |
! Input/output: |
166 |
! A(N,N) : Matrix to be diagonalized, store in a (Nreal,Nreal) matrix |
167 |
! Destroyed in output. |
168 |
! Output: |
169 |
! V(N,N) : Eigenvectors, stored in V(NReal, NReal) |
170 |
! D(N) : Eigenvalues, stored in D(NReal) |
171 |
! |
172 |
|
173 |
Use Vartypes |
174 |
|
175 |
IMPLICIT NONE |
176 |
INTEGER(KINT), parameter :: max_it=500 |
177 |
REAL(KREAL), ALLOCATABLE :: B(:),Z(:) |
178 |
|
179 |
INTEGER(KINT) :: N,NReal |
180 |
REAL(KREAL) :: A(NReal,NReal) |
181 |
REAL(KREAL) :: V(Nreal,Nreal),D(Nreal) |
182 |
|
183 |
INTEGER(KINT) :: I, J,IP, IQ, NROT |
184 |
REAL(KREAL) :: SM, H, Tresh, G, T, Theta, C, S, Tau |
185 |
|
186 |
allocate(B(N),Z(N)) |
187 |
|
188 |
DO IP=1,N |
189 |
DO IQ=1,N |
190 |
V(IP,IQ)=0. |
191 |
END DO |
192 |
V(IP,IP)=1. |
193 |
END DO |
194 |
DO IP=1,N |
195 |
B(IP)=A(IP,IP) |
196 |
D(IP)=B(IP) |
197 |
Z(IP)=0. |
198 |
END DO |
199 |
NROT=0 |
200 |
DO I=1,max_it |
201 |
SM=0. |
202 |
DO IP=1,N-1 |
203 |
DO IQ=IP+1,N |
204 |
SM=SM+ABS(A(IP,IQ)) |
205 |
END DO |
206 |
END DO |
207 |
IF(SM.EQ.0.) GOTO 100 |
208 |
IF(I.LT.4)THEN |
209 |
TRESH=0.2*SM/N**2 |
210 |
ELSE |
211 |
TRESH=0. |
212 |
ENDIF |
213 |
DO IP=1,N-1 |
214 |
DO IQ=IP+1,N |
215 |
G=100.*ABS(A(IP,IQ)) |
216 |
IF((I.GT.4).AND.(ABS(D(IP))+G.EQ.ABS(D(IP))) & |
217 |
.AND.(ABS(D(IQ))+G.EQ.ABS(D(IQ))))THEN |
218 |
A(IP,IQ)=0. |
219 |
ELSE IF(ABS(A(IP,IQ)).GT.TRESH)THEN |
220 |
H=D(IQ)-D(IP) |
221 |
IF(ABS(H)+G.EQ.ABS(H))THEN |
222 |
T=A(IP,IQ)/H |
223 |
ELSE |
224 |
THETA=0.5*H/A(IP,IQ) |
225 |
T=1./(ABS(THETA)+SQRT(1.+THETA**2)) |
226 |
IF(THETA.LT.0.)T=-T |
227 |
ENDIF |
228 |
C=1./SQRT(1+T**2) |
229 |
S=T*C |
230 |
TAU=S/(1.+C) |
231 |
H=T*A(IP,IQ) |
232 |
Z(IP)=Z(IP)-H |
233 |
Z(IQ)=Z(IQ)+H |
234 |
D(IP)=D(IP)-H |
235 |
D(IQ)=D(IQ)+H |
236 |
A(IP,IQ)=0. |
237 |
DO J=1,IP-1 |
238 |
G=A(J,IP) |
239 |
H=A(J,IQ) |
240 |
A(J,IP)=G-S*(H+G*TAU) |
241 |
A(J,IQ)=H+S*(G-H*TAU) |
242 |
END DO |
243 |
DO J=IP+1,IQ-1 |
244 |
G=A(IP,J) |
245 |
H=A(J,IQ) |
246 |
A(IP,J)=G-S*(H+G*TAU) |
247 |
A(J,IQ)=H+S*(G-H*TAU) |
248 |
END DO |
249 |
DO J=IQ+1,N |
250 |
G=A(IP,J) |
251 |
H=A(IQ,J) |
252 |
A(IP,J)=G-S*(H+G*TAU) |
253 |
A(IQ,J)=H+S*(G-H*TAU) |
254 |
END DO |
255 |
DO J=1,N |
256 |
G=V(J,IP) |
257 |
H=V(J,IQ) |
258 |
V(J,IP)=G-S*(H+G*TAU) |
259 |
V(J,IQ)=H+S*(G-H*TAU) |
260 |
END DO |
261 |
NROT=NROT+1 |
262 |
ENDIF |
263 |
END DO |
264 |
END DO |
265 |
DO IP=1,N |
266 |
B(IP)=B(IP)+Z(IP) |
267 |
D(IP)=B(IP) |
268 |
Z(IP)=0. |
269 |
END DO |
270 |
END DO |
271 |
write(6,*) max_it,' iterations should never happen' |
272 |
STOP |
273 |
100 DEALLOCATE(B,Z) |
274 |
RETURN |
275 |
END SUBROUTINE JACOBI |
276 |
! |
277 |
!============================================================ |
278 |
!c |
279 |
!c ++ Sort Eigenvectors in ascending eigenvalues order |
280 |
!c |
281 |
!c============================================================ |
282 |
!c |
283 |
SUBROUTINE trie(nb_niv,ene,psi,nreal) |
284 |
|
285 |
! |
286 |
! Input: |
287 |
! Nb_niv : Dimension of Psi |
288 |
! NReal : Actual dimensions of Psi, Ene |
289 |
! Input/output: |
290 |
! Psi(N,N) : Eigvenvectors, changed in output. Stored in a (Nreal,Nreal) matrix |
291 |
! Ene(N) : Eigenvalues, changed in output. Stored in Ene(NReal) |
292 |
! |
293 |
!!!!!!!!!!!!!!!! |
294 |
|
295 |
Use VarTypes |
296 |
IMPLICIT NONE |
297 |
|
298 |
integer(KINT) :: i, j, k, nb_niv, nreal |
299 |
real(KREAL) :: ene(nreal),psi(nreal,nreal) |
300 |
real(KREAL) :: a |
301 |
|
302 |
|
303 |
!!!!!!!!!!!!!!!! |
304 |
|
305 |
|
306 |
DO i=1,nb_niv |
307 |
DO j=i+1,nb_niv |
308 |
IF (ene(i) .GT. ene(j)) THEN |
309 |
! permutation |
310 |
a=ene(i) |
311 |
ene(i)=ene(j) |
312 |
ene(j)=a |
313 |
|
314 |
DO k=1,nb_niv |
315 |
a=psi(k,i) |
316 |
psi(k,i)=psi(k,j) |
317 |
psi(k,j)=a |
318 |
END DO |
319 |
END IF |
320 |
END DO |
321 |
END DO |
322 |
|
323 |
END SUBROUTINE trie |
324 |
|