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
|
#include "log-fourier.hpp"
LogarithmicFourierTransform::LogarithmicFourierTransform(unsigned N, Real k, Real Δτ, unsigned pad) : N(N), pad(pad), k(k), Δτ(Δτ) {
τₛ = -0.5 * N;
ωₛ = -0.5 * N;
sₛ = -0.5 * pad * N;
a = reinterpret_cast<Complex*>(fftw_alloc_complex(pad*N));
â = reinterpret_cast<Complex*>(fftw_alloc_complex(pad*N));
fftw_import_wisdom_from_filename("fftw.wisdom");
a_to_â = fftw_plan_dft_1d(pad*N, reinterpret_cast<fftw_complex*>(a), reinterpret_cast<fftw_complex*>(â), FFTW_BACKWARD, 0);
â_to_a = fftw_plan_dft_1d(pad*N, reinterpret_cast<fftw_complex*>(â), reinterpret_cast<fftw_complex*>(a), FFTW_BACKWARD, 0);
fftw_export_wisdom_to_filename("fftw.wisdom");
}
LogarithmicFourierTransform::~LogarithmicFourierTransform() {
fftw_destroy_plan(a_to_â);
fftw_destroy_plan(â_to_a);
fftw_free(a);
fftw_free(â);
fftw_cleanup();
}
Real LogarithmicFourierTransform::τ(unsigned n) const {
return Δτ * (n + τₛ);
}
Real LogarithmicFourierTransform::ω(unsigned n) const {
return Δτ * (n + ωₛ);
}
Real LogarithmicFourierTransform::s(unsigned n) const {
return (n + sₛ) * 2*M_PI / (pad * N * Δτ);
}
Real LogarithmicFourierTransform::t(unsigned n) const {
return exp(τ(n));
}
Real LogarithmicFourierTransform::ν(unsigned n) const {
return exp(ω(n));
}
Complex Γ(Complex z) {
gsl_sf_result logΓ;
gsl_sf_result argΓ;
gsl_sf_lngamma_complex_e(z.real(), z.imag(), &logΓ, &argΓ);
return exp(logΓ.val + 1i * argΓ.val);
}
std::vector<Complex> LogarithmicFourierTransform::fourier(const std::vector<Real>& c, bool symmetric) {
std::vector<Complex> ĉ(N);
std::vector<Real> σs = {1};
/* c is either even or zero for negative arguments */
if (symmetric){
σs.push_back(-1);
}
for (Real σ : σs) {
for (unsigned n = 0; n < pad*N; n++) {
if (n < N) {
a[n] = c[n] * exp((1 - k) * τ(n));
} else {
a[n] = 0;
}
}
fftw_execute(a_to_â);
for (unsigned n = 0; n < pad*N; n++) {
â[(pad*N / 2 + n) % (pad*N)] *= pow(1i * σ, 1i * s(n) - k) * Γ(k - 1i * s(n));
}
fftw_execute(â_to_a);
for (unsigned n = 0; n < N; n++) {
ĉ[n] += exp(-k * ω(n)) * a[(pad - 1)*N+n] / (Real)(pad*N);
}
}
return ĉ;
}
std::vector<Real> LogarithmicFourierTransform::inverse(const std::vector<Complex>& ĉ) {
std::vector<Real> c(N);
std::vector<Real> σs = {1, -1};
for (Real σ : σs) {
for (unsigned n = 0; n < pad * N; n++) {
if (n < N) {
a[n] = ĉ[n] * exp((1 - k) * ω(n));
} else {
a[n] = 0;
}
}
fftw_execute(a_to_â);
for (unsigned n = 0; n < pad*N; n++) {
â[(pad*N / 2 + n) % (pad*N)] *= pow(-1i * σ, 1i * s(n) - k) * Γ(k - 1i * s(n));
}
fftw_execute(â_to_a);
for (unsigned n = 0; n < N; n++) {
c[n] += exp(-k * τ(n)) * a[(pad - 1)*N+n].real() / (Real)(pad*N) / (2 * M_PI);
}
}
return c;
}
|