-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathapfft.f90
324 lines (251 loc) · 8.51 KB
/
apfft.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
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
!+
! module apfft_mod
!
! This module implements the All Phase FFT method for obtaining accurate phase from signal data.
!-
module apfft_mod
use physical_constants
implicit none
contains
!--------------------------------------------------------------------------------------------------
!--------------------------------------------------------------------------------------------------
!--------------------------------------------------------------------------------------------------
!+
! subroutine apfft_corr(rdata_in, bounds, window, phase, amp, freq, diag)
!
! For real signal rdata_in, computes phase, frequency, and amplitude
! of peak found within bounds. Algorithm is corrected all-phase FFT and should.
!
! This routine finds only one peak: the largest amplitude within the bound. Signals with multiple
! components can be investigated by varying bounds appropriately.
!
! Input:
! rdata_in(:) -- real(rp): signal data.
! bounds(1:2) -- real(rp): range within which to search for peak.
! window -- character(3): 'rec' or 'han' for rectangular or Hann window.
! diag -- integer, optional: causes low-level routine apfft_ext to produce
! a fort.X file where X=9000+fid containing diag information.
! Output:
! phase -- real(rp): phase of peak found in signal.
! freq -- real(rp): frequency of peak
! amp -- real(rp): amplitude of peak
!-
subroutine apfft_corr(rdata_in, bounds, window, phase, amp, freq, diag)
use, intrinsic :: ieee_arithmetic
implicit none
real(rp) rdata_in(:)
real(rp) phase, amp, freq
real(rp), optional :: bounds(2)
character(3) window
integer, optional :: diag
real(rp), allocatable :: rdata(:)
real(rp), allocatable :: u1(:), u2(:)
real(rp) phase_u1, amp_u1, freq_u1
real(rp) phase_u2, amp_u2, freq_u2
real(rp) d
integer i,N, ndata
character(*), parameter :: r_name = 'apfft_corr'
ndata = size(rdata_in)
N = floor((ndata+1.0)/3)
ndata = 3*N-1
allocate(rdata(ndata))
!remove zero offset
rdata = rdata_in(1:ndata) - sum(rdata_in(1:ndata))/size(rdata_in(1:ndata))
allocate(u1(2*N-1))
allocate(u2(2*N-1))
u1 = rdata(1:2*N-1)
u2 = rdata(N+1:3*N-1)
call apfft_ext(u1, bounds, window, phase_u1, amp_u1, freq_u1, diag)
call apfft_ext(u2, bounds, window, phase_u2, amp_u2, freq_u2, diag)
d = (phase_u2-phase_u1)/twopi
if(d .gt. 0.5) then
d = d - 1.0d0
elseif(d .le. -0.5) then
d = d + 1.0d0
endif
freq = freq_u1 + (1.0d0*d)/N
phase = 2.0d0*phase_u1 - phase_u2
if(phase .lt. 0) then
phase = phase + twopi
elseif(phase .gt. twopi) then
phase = phase - twopi
endif
if(window == 'rec') then
if( ieee_is_finite(1.0d0/sin(pi*d)) ) then
amp = (pi*d/sin(pi*d))**2 * amp_u1 * 2.0
else
amp = 2.0d0 * amp_u1
endif
elseif(window == 'han') then
if( ieee_is_finite(1.0d0/sin(pi*d)) ) then
amp = (pi*d*(1.0d0-d*d)/sin(pi*d))**2 * amp_u1 * 2.0
else
amp = 2.0d0 * amp_u1
endif
endif
deallocate(rdata)
deallocate(u1)
deallocate(u2)
end subroutine apfft_corr
!--------------------------------------------------------------------------------------------------
!--------------------------------------------------------------------------------------------------
!--------------------------------------------------------------------------------------------------
!+
! subroutine apfft(rdata_in, bounds, window, phase, diag)
!
! Implements the All Phase FFT method for obtaining accurate phase from signal data.
!
! The signal data is truncated to an odd length, and the phase is relative to the central point.
!
! Input:
! rdata_in(:) - real(rp): signal data.
! bounds(2) - real(rp): Fractional. Boundaries within which peak is searched.
! window - character(3): 'rec' or 'han'
! diag - integer, optional :: If present and positive, produced diag file containing spectrum.
!
! Output:
!
! phase - real(rp), Phase of found peak, as determined by ApFFT.
!-
subroutine apfft(rdata_in, bounds, window, phase, diag)
use fgsl
implicit none
real(rp) rdata_in(:)
real(rp) bounds(2)
character(3) window
real(rp) phase
real(rp) amp, freq
real(rp) rdata(size(rdata_in))
integer, optional :: diag
rdata = rdata_in - 1.0d0*sum(rdata_in)/size(rdata_in)
call apfft_ext(rdata, bounds, window, phase, amp, freq, diag)
end subroutine apfft
!--------------------------------------------------------------------------------------------------
!--------------------------------------------------------------------------------------------------
!--------------------------------------------------------------------------------------------------
!+
! subroutine apfft_ext(rdata,bounds, window, phase, amp, freq, diag)
!
! Implements the All Phase FFT method for obtaining accurate phase from signal data.
!
! This "extended" apfft subroutine returns the amplitudes and frequency as well, for use
! by the corrected apfft subroutine in this module.
!
! Input:
! rdata_in(:) - real(rp): signal data.
!
! Output:
! phase - real(rp), Phase dominant peak, as determined by ApFFT.
!-
subroutine apfft_ext(rdata, bounds, window, phase, amp, freq, diag)
use fgsl
implicit none
real(rp) rdata(:)
real(rp) bounds(2)
character(3) window
real(rp) phase, amp, freq
integer, optional :: diag
real(rp) wc(size(rdata))
integer fid
integer ndata, N, max_ix
integer i, j, lb, ub
real(rp) alpha, beta, gamma, p
complex(rp), allocatable :: c_Xap(:)
type(fgsl_fft_complex_wavetable) :: wavetable
type(fgsl_fft_complex_workspace) :: work
integer(fgsl_int) :: status
real(rp), allocatable :: apwindow(:)
character(4) fid_str
character(*), parameter :: r_name = 'apfft_ext'
ndata = size(rdata)
N = floor((ndata+1.0d0)/2.0d0)
allocate(c_Xap(N))
!Produce the all-phase vector from signal data.
if(window == 'rec') then
!rec-rec window
c_Xap(1) = cmplx(N*rdata(N),0.0d0)
do i=2,N
c_Xap(i) = cmplx((N-i+1)*rdata(N+i-1) + (i-1)*rdata(i-1),0.0d0)
enddo
c_Xap = c_Xap / N / N
elseif( window == 'han') then
allocate(apwindow(ndata))
call hanhan(N,apwindow)
rdata = apwindow*rdata
deallocate(apwindow)
c_Xap(1) = cmplx(rdata(N),0.0d0)
do i=2,N
c_Xap(i) = cmplx(rdata(N+i-1) + rdata(i-1),0.0d0)
enddo
c_Xap = c_Xap
else
write(*,*) "Unknown window passed to apfft_ext: ", window
call err_exit()
endif
wavetable = fgsl_fft_complex_wavetable_alloc(int(N,fgsl_size_t))
work = fgsl_fft_complex_workspace_alloc(int(N,fgsl_size_t))
status = fgsl_fft_complex_backward(c_Xap, 1_fgsl_size_t, int(N,fgsl_size_t), wavetable, work)
call fgsl_fft_complex_workspace_free(work)
call fgsl_fft_complex_wavetable_free(wavetable)
if(present(diag)) then
if( diag .gt. 0) then
fid = 9000 + diag
write(fid_str,'(i4)') fid
open(fid, file = 'apfft_'//fid_str//'.diag')
write(fid,*) '#apfft_ext diag file'
write(fid,*) '#freq abs angle'
do i=1,N
write(fid,'(3f15.6)') (i-1.0d0)/N, abs(c_Xap(i)), twopi_atan2(aimag(c_Xap(i)),real(c_Xap(i)))
enddo
write(fid,*)
write(fid,*)
close(fid)
endif
endif
lb = nint(N*bounds(1)) + 1
ub = nint(N*bounds(2)) + 1
max_ix = maxloc(abs(c_Xap(lb:ub)),1) + lb - 1
phase = twopi_atan2(aimag(c_Xap(max_ix)),real(c_Xap(max_ix)))
amp = abs(c_Xap(max_ix))
freq = (max_ix-1.0d0)/N
deallocate(c_Xap)
!----------------------------------------------------------------------
contains
function twopi_atan2(Y,X) result(arg)
use physical_constants
real(rp) Y, X, arg
arg = atan2(Y,X)
if(arg < 0) arg = arg+twopi
arg = twopi - arg
end function twopi_atan2
end subroutine apfft_ext
!--------------------------------------------------------------------------------------------------
!--------------------------------------------------------------------------------------------------
!--------------------------------------------------------------------------------------------------
subroutine hanhan(N, hh)
integer N
real(rp) hh(:) !length must be 2*N-1
real(rp) h(N)
integer i
!
h = han(N)
!convolution of two Hann windows
do i=1,N
hh(i) = sum( h(N-i+1:N) * h(1:i) )
hh(2*N-i) = hh(i)
enddo
hh = hh
end subroutine hanhan
!--------------------------------------------------------------------------------------------------
!--------------------------------------------------------------------------------------------------
!--------------------------------------------------------------------------------------------------
function han(N)
implicit none
integer i, N
real(rp) han(N)
!
do i=1,N
han(i) = sin(pi*i/(N-1.0d0))**2 / (N-1.0d0) * 2.0d0
enddo
end function han
end module