From d93dc5dba5d9e00298038c3dac3b51aa554d2ca8 Mon Sep 17 00:00:00 2001
From: Jaron Kent-Dobias <jaron@kent-dobias.com>
Date: Sat, 10 May 2025 13:46:24 -0300
Subject: New time scaling

---
 log-fourier.cpp | 29 +++++++++++++++--------------
 1 file changed, 15 insertions(+), 14 deletions(-)

(limited to 'log-fourier.cpp')

diff --git a/log-fourier.cpp b/log-fourier.cpp
index e60439a..d5ac17d 100644
--- a/log-fourier.cpp
+++ b/log-fourier.cpp
@@ -4,9 +4,9 @@
 #include <fstream>
 #include <types.hpp>
 
-LogarithmicFourierTransform::LogarithmicFourierTransform(unsigned N, Real k, Real Δτ, unsigned pad, Real shift) : N(N), pad(pad), k(k), Δτ(Δτ) {
-  τₛ = -shift * N;
-  ωₛ = -(1-shift) * N;
+LogarithmicFourierTransform::LogarithmicFourierTransform(unsigned N, Real k, Real Δτ, unsigned pad, Real shift) : N(N), pad(pad), k(k), Δτ(Δτ), shift(shift) {
+  τₛ = -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));
@@ -37,11 +37,11 @@ Real LogarithmicFourierTransform::s(unsigned n) const {
 }
 
 Real LogarithmicFourierTransform::t(unsigned n) const {
-  return exp(τ(n));
+  return std::exp(τ(n)) / shift;
 }
 
 Real LogarithmicFourierTransform::ν(unsigned n) const {
-  return exp(ω(n));
+  return std::exp(ω(n)) * shift;
 }
 
 Complex Γ(Complex z) {
@@ -50,7 +50,7 @@ Complex Γ(Complex z) {
 
   gsl_sf_lngamma_complex_e((double)z.real(), (double)z.imag(), &logΓ, &argΓ);
 
-  return exp((Real)logΓ.val + II * (Real)argΓ.val);
+  return std::exp((Real)logΓ.val + II * (Real)argΓ.val);
 }
 
 std::vector<Complex> LogarithmicFourierTransform::fourier(const std::vector<Real>& c, bool symmetric) {
@@ -63,16 +63,16 @@ std::vector<Complex> LogarithmicFourierTransform::fourier(const std::vector<Real
   for (Real σ : σs) {
     for (unsigned n = 0; n < pad*N; n++) {
       if (n < N) {
-        a[n] = c[n] * exp((1 - k) * τ(n));
+        a[n] = c[n] * std::exp((1 - k) * τ(n));
       } else if (n >= (pad - 1) * N) {
-        a[n] = c[pad*N-n-1] * exp((1 - k) * τ(pad*N-n-1));
+        a[n] = c[pad*N-n-1] * std::exp((1 - k) * τ(pad*N-n-1));
       } else {
         a[n] = 0;
       }
     }
     FFTW_EXECUTE(a_to_â);
     for (unsigned n = 0; n < pad*N; n++) {
-      â[(pad*N / 2 + n) % (pad*N)] *= std::exp(II*(0.5 * N + τₛ) * s(n) / Δτ) * std::pow(II * σ, II * s(n) - k) * Γ(k - II * s(n));
+      â[(pad*N / 2 + n) % (pad*N)] *= std::pow(II * σ, II * s(n) - k) * Γ(k - II * s(n));
     }
     FFTW_EXECUTE(â_to_a);
     for (unsigned n = 0; n < N; n++) {
@@ -82,6 +82,8 @@ std::vector<Complex> LogarithmicFourierTransform::fourier(const std::vector<Real
 
   for (unsigned n = 0; n < N; n++) {
     ĉ[n] -= ĉ[N - 1];
+    if (symmetric) ĉ[n] = ĉ[n].real();
+    ĉ[n] /= shift;
   }
 
   return ĉ;
@@ -95,14 +97,14 @@ std::vector<Real> LogarithmicFourierTransform::inverse(const std::vector<Complex
       if (n < N) {
         a[n] = (ĉ[n].real() + II * σ * ĉ[n].imag()) * std::exp((1 - k) * ω(n));
       } else if (n >= (pad - 1) * N) {
-        a[n] = (ĉ[pad*N-n-1].real() + II * σ * ĉ[pad*N-n-1].imag()) * std::exp((1 - k) * ω(pad*N-n-1));
+        a[n] = (ĉ[pad*N-n-1].real() - II * σ * ĉ[pad*N-n-1].imag()) * std::exp((1 - k) * ω(pad*N-n-1));
       } else {
         a[n] = 0;
       }
     }
     FFTW_EXECUTE(a_to_â);
     for (unsigned n = 0; n < pad*N; n++) {
-      â[(pad*N / 2 + n) % (pad*N)] *= std::exp(-II*(0.5 * N + τₛ) * s(n) / Δτ) * std::pow(-II * σ, II * s(n) - k) * Γ(k - II * s(n));
+      â[(pad*N / 2 + n) % (pad*N)] *= std::pow(-II * σ, II * s(n) - k) * Γ(k - II * s(n));
     }
     FFTW_EXECUTE(â_to_a);
     for (unsigned n = 0; n < N; n++) {
@@ -112,6 +114,7 @@ std::vector<Real> LogarithmicFourierTransform::inverse(const std::vector<Complex
 
   for (unsigned n = 0; n < N; n++) {
     c[n] -= c[N - 1];
+    c[n] *= shift;
   }
 
   return c;
@@ -189,11 +192,9 @@ Real estimateZ(LogarithmicFourierTransform& fft, const std::vector<Real>& C, con
 
 Real energy(const LogarithmicFourierTransform& fft, std::vector<Real>&  C, const std::vector<Real>& R, unsigned p, unsigned s, Real λ, Real β) {
   unsigned n₀ = 0;
-  /*
   for (unsigned n = 0; n < C.size(); n++) {
     if (C[n] > 1 || R[n] > 1) n₀ = n % 2 == 0 ? n / 2 : (n + 1) / 2;
   }
-  */
   Real E = fft.t(2*n₀) * df(λ, p, s, 1);
   for (unsigned n = n₀; n < C.size()/2-1; n++) {
     Real R₂ₙ   = R[2*n];
@@ -203,7 +204,7 @@ Real energy(const LogarithmicFourierTransform& fft, std::vector<Real>&  C, const
     Real C₂ₙ₊₁ = C[2*n+1];
     Real C₂ₙ₊₂ = C[2*n+2];
 
-    //if (C₂ₙ₊₂ < 0 || R₂ₙ₊₂ < 0) break;
+    if (C₂ₙ₊₂ < 0 || R₂ₙ₊₂ < 0) break;
 
     Real h₂ₙ   = fft.t(2*n+1) - fft.t(2*n);
     Real h₂ₙ₊₁ = fft.t(2*n+2) - fft.t(2*n+1);
-- 
cgit v1.2.3-70-g09d2