Commit b85eed0171def3a9c57818323d3494afee5fc0ed
1 parent
19e954a45e
Exists in
master
and in
3 other branches
Added sorting of eigenvalues and vectors using an insertion sort algo to better
match imsl behavior git-svn-id: https://lxsd.femto-st.fr/svn/fvn@52 b657c933-2333-4658-acf2-d3c7c2708721
Showing 1 changed file with 73 additions and 0 deletions Side-by-side Diff
fvn_linear/fvn_linear.f90
| ... | ... | @@ -1025,6 +1025,67 @@ |
| 1025 | 1025 | ! |
| 1026 | 1026 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
| 1027 | 1027 | |
| 1028 | +! August 2009 | |
| 1029 | +! William Daniau | |
| 1030 | +! Adding sorting of eigenvalues and vectors | |
| 1031 | +subroutine fvn_z_sort_eigen(d,v,vec) | |
| 1032 | +! this routine takes in input : | |
| 1033 | +! v : a vector containing unsorted eigenvalues | |
| 1034 | +! vec : a matrix where vec(:,j) is the eigenvector corresponding to v(j) | |
| 1035 | +! | |
| 1036 | +! At the end of subroutine the eigenvalues are sorted by decreasing order of | |
| 1037 | +! modulus so as the vectors. | |
| 1038 | +implicit none | |
| 1039 | +integer :: d,i,j | |
| 1040 | +complex(8),dimension(d) :: v | |
| 1041 | +complex(8),dimension(d,d) :: vec | |
| 1042 | +complex(8) :: cur_v | |
| 1043 | +complex(8), dimension(d) :: cur_vec | |
| 1044 | + | |
| 1045 | +do i=2,d | |
| 1046 | + cur_v=v(i) | |
| 1047 | + cur_vec=vec(:,i) | |
| 1048 | + j=i-1 | |
| 1049 | + do while( (j>=1) .and. (abs(v(j)) < abs(cur_v)) ) | |
| 1050 | + v(j+1)=v(j) | |
| 1051 | + vec(:,j+1)=vec(:,j) | |
| 1052 | + j=j-1 | |
| 1053 | + end do | |
| 1054 | + v(j+1)=cur_v | |
| 1055 | + vec(:,j+1)=cur_vec | |
| 1056 | +end do | |
| 1057 | + | |
| 1058 | +end subroutine | |
| 1059 | + | |
| 1060 | +subroutine fvn_c_sort_eigen(d,v,vec) | |
| 1061 | +! this routine takes in input : | |
| 1062 | +! v : a vector containing unsorted eigenvalues | |
| 1063 | +! vec : a matrix where vec(:,j) is the eigenvector corresponding to v(j) | |
| 1064 | +! | |
| 1065 | +! At the end of subroutine the eigenvalues are sorted by decreasing order of | |
| 1066 | +! modulus so as the vectors. | |
| 1067 | +implicit none | |
| 1068 | +integer :: d,i,j | |
| 1069 | +complex(4),dimension(d) :: v | |
| 1070 | +complex(4),dimension(d,d) :: vec | |
| 1071 | +complex(4) :: cur_v | |
| 1072 | +complex(4), dimension(d) :: cur_vec | |
| 1073 | + | |
| 1074 | +do i=2,d | |
| 1075 | + cur_v=v(i) | |
| 1076 | + cur_vec=vec(:,i) | |
| 1077 | + j=i-1 | |
| 1078 | + do while( (j>=1) .and. (abs(v(j)) < abs(cur_v)) ) | |
| 1079 | + v(j+1)=v(j) | |
| 1080 | + vec(:,j+1)=vec(:,j) | |
| 1081 | + j=j-1 | |
| 1082 | + end do | |
| 1083 | + v(j+1)=cur_v | |
| 1084 | + vec(:,j+1)=cur_vec | |
| 1085 | +end do | |
| 1086 | + | |
| 1087 | +end subroutine | |
| 1088 | + | |
| 1028 | 1089 | subroutine fvn_s_matev(d,a,evala,eveca,status) |
| 1029 | 1090 | ! |
| 1030 | 1091 | ! integer d (in) : matrice rank |
| ... | ... | @@ -1098,6 +1159,9 @@ |
| 1098 | 1159 | deallocate(wr) |
| 1099 | 1160 | deallocate(wc_a) |
| 1100 | 1161 | |
| 1162 | + ! sorting | |
| 1163 | + call fvn_c_sort_eigen(d,evala,eveca) | |
| 1164 | + | |
| 1101 | 1165 | end subroutine |
| 1102 | 1166 | |
| 1103 | 1167 | subroutine fvn_d_matev(d,a,evala,eveca,status) |
| ... | ... | @@ -1174,6 +1238,9 @@ |
| 1174 | 1238 | deallocate(wr) |
| 1175 | 1239 | deallocate(wc_a) |
| 1176 | 1240 | |
| 1241 | + ! sorting | |
| 1242 | + call fvn_z_sort_eigen(d,evala,eveca) | |
| 1243 | + | |
| 1177 | 1244 | end subroutine |
| 1178 | 1245 | |
| 1179 | 1246 | subroutine fvn_c_matev(d,a,evala,eveca,status) |
| ... | ... | @@ -1222,6 +1289,9 @@ |
| 1222 | 1289 | deallocate(work) |
| 1223 | 1290 | deallocate(wc_a) |
| 1224 | 1291 | |
| 1292 | + ! sorting | |
| 1293 | + call fvn_c_sort_eigen(d,evala,eveca) | |
| 1294 | + | |
| 1225 | 1295 | end subroutine |
| 1226 | 1296 | |
| 1227 | 1297 | subroutine fvn_z_matev(d,a,evala,eveca,status) |
| ... | ... | @@ -1268,6 +1338,9 @@ |
| 1268 | 1338 | deallocate(rwork) |
| 1269 | 1339 | deallocate(work) |
| 1270 | 1340 | deallocate(wc_a) |
| 1341 | + | |
| 1342 | + ! sorting | |
| 1343 | + call fvn_z_sort_eigen(d,evala,eveca) | |
| 1271 | 1344 | |
| 1272 | 1345 | end subroutine |
| 1273 | 1346 |