-
Notifications
You must be signed in to change notification settings - Fork 0
/
correlation_integral.f90
152 lines (127 loc) · 3.62 KB
/
correlation_integral.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
! Correlation Integral caluclations done fast with openmp
! Copyright (C) 2021 Victor Azizi - [email protected]
! This program is free software: you can redistribute it and/or modify
! it under the terms of the GNU General Public License as published by
! the Free Software Foundation, either version 3 of the License, or
! (at your option) any later version.
!
! This program 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 General Public License for more details.
!
! You should have received a copy of the GNU General Public License
! along with this program. If not, see <https://www.gnu.org/licenses/>.
subroutine manhattan(data, dims, r, cd)
use iso_fortran_env
implicit none
real(real32) , intent(in) :: data(:)
integer , intent(in) :: dims
real(real32) , intent(in) :: r(:)
real(real32) , intent(out):: cd(size(r))
integer(int64) :: cd_c(size(r))
integer(int64), allocatable:: cd_p(:,:)
integer :: L, L_2, i, j, k, datasize, rsize
real(real32) :: dist
cd_c = 0
datasize = size(data)
rsize = size(r)
L = datasize - dims
L_2 = L/2
allocate(cd_p(rsize,L))
cd_p = 0
!$OMP PARALLEL DO PRIVATE(dist) schedule(guided)
do i=L,1,-1
do j=i+1,L
dist = 0
do k=1,dims
dist = dist + abs ( data(i+k) - data(j+k) )
enddo
do k=1,rsize
if (r(k) > dist) then
cd_p(k,i) = cd_p(k,i) + 1
endif
enddo
enddo
enddo
!$OMP END PARALLEL DO
cd_c = sum( cd_p,2)
cd = real(real(cd_c, real64) / (real(L, real64) * (real(L, real64) - 1. ) / 2.),real32)
deallocate(cd_p)
end subroutine
subroutine euclidean(data, dims, r, cd)
use iso_fortran_env
implicit none
real(real32) , intent(in) :: data(:)
integer , intent(in) :: dims
real(real32) , intent(in) :: r(:)
real(real32) , intent(out):: cd(size(r))
integer(int64) :: cd_c(size(r))
integer(int64), allocatable:: cd_p(:,:)
integer :: L, L_2, i, j, k, datasize, rsize
real(real32) :: dist
real(real32) :: r_p2(size(r))
cd_c = 0
datasize = size(data)
rsize = size(r)
L = datasize - dims
L_2 = L/2
r_p2 = r**2
allocate(cd_p(rsize,L))
cd_p = 0
!$OMP PARALLEL DO PRIVATE(dist) schedule(guided)
do i=L,1,-1
do j=i+1,L
dist = 0
do k=1,dims
dist = dist + (data(i+k) - data(j+k))**2
enddo
do k=1,rsize
if (r_p2(k) > dist) then
cd_p(k,i) = cd_p(k,i) + 1
endif
enddo
enddo
enddo
!$OMP END PARALLEL DO
cd_c = sum( cd_p,2)
cd = real(real(cd_c, real64) / (real(L, real64) * (real(L, real64) - 1. ) / 2.),real32)
deallocate(cd_p)
end subroutine
subroutine chebyshev(data, dims, r, cd)
use iso_fortran_env
implicit none
real(real32) , intent(in) :: data(:)
integer , intent(in) :: dims
real(real32) , intent(in) :: r(:)
real(real32) , intent(out):: cd(size(r))
integer(int64) :: cd_c(size(r))
integer(int64), allocatable:: cd_p(:,:)
integer :: L, L_2, i, j, k, datasize, rsize
real(real32) :: dist
cd_c = 0
datasize = size(data)
rsize = size(r)
L = datasize - dims
L_2 = L/2
allocate(cd_p(rsize,L))
cd_p = 0
!$OMP PARALLEL DO PRIVATE(dist) schedule(guided)
do i=L,1,-1
do j=i+1,L
dist = 0
do k=1,dims
dist = max(dist, abs(data(i+k) - data(j+k)))
enddo
do k=1,rsize
if (r(k) > dist) then
cd_p(k,i) = cd_p(k,i) + 1
endif
enddo
enddo
enddo
!$OMP END PARALLEL DO
cd_c = sum( cd_p,2)
cd = real(real(cd_c, real64) / (real(L, real64) * (real(L, real64) - 1. ) / 2.),real32)
deallocate(cd_p)
end subroutine