#include "log-fourier.hpp"
#include "p-spin.hpp"
#include <complex>
#include <fstream>

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) {
        if (σ < 0) {
          a[n] = std::conj(ĉ[n]) * exp((1 - k) * ω(n));
        } else {
          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;
}

std::string logFourierFile(std::string prefix, unsigned p, unsigned s, Real λ, Real τ₀, Real β, unsigned log2n, Real Δτ, Real k) {
  return prefix + "_" + std::to_string(p) + "_" + std::to_string(s) + "_" + std::to_string(λ) + "_" + std::to_string(τ₀) + "_" + std::to_string(β) + "_" + std::to_string(log2n) + "_" + std::to_string(Δτ)  + "_" + std::to_string(k) + ".dat";
}

void logFourierSave(const std::vector<Real>& C, const std::vector<Real>& R, const std::vector<Complex>& Ct, const std::vector<Complex>& Rt, unsigned p, unsigned s, Real λ, Real τ₀, Real β, unsigned log2n, Real Δτ, Real k) {
    unsigned N = pow(2, log2n);
    std::ofstream outfile(logFourierFile("C", p, s, λ, τ₀, β, log2n, Δτ, k), std::ios::out | std::ios::binary);
    outfile.write((const char*)(C.data()), N * sizeof(Real));
    outfile.close();

    std::ofstream outfileCt(logFourierFile("Ct", p, s, λ, τ₀, β, log2n, Δτ, k), std::ios::out | std::ios::binary);
    outfileCt.write((const char*)(Ct.data()), N * sizeof(Complex));
    outfileCt.close();

    std::ofstream outfileR(logFourierFile("R", p, s, λ, τ₀, β, log2n, Δτ, k), std::ios::out | std::ios::binary);
    outfileR.write((const char*)(R.data()), N * sizeof(Real));
    outfileR.close();

    std::ofstream outfileRt(logFourierFile("Rt", p, s, λ, τ₀, β, log2n, Δτ, k), std::ios::out | std::ios::binary);
    outfileRt.write((const char*)(Rt.data()), N * sizeof(Complex));
    outfileRt.close();
}

bool logFourierLoad(std::vector<Real>& C, std::vector<Real>& R, std::vector<Complex>& Ct, std::vector<Complex>& Rt, unsigned p, unsigned s, Real λ, Real τ₀, Real β, unsigned log2n, Real Δτ, Real k) {
  std::ifstream cfile(logFourierFile("C", p, s, λ, τ₀, β, log2n, Δτ, k), std::ios::binary);
  std::ifstream rfile(logFourierFile("R", p, s, λ, τ₀, β, log2n, Δτ, k), std::ios::binary);
  std::ifstream ctfile(logFourierFile("Ct", p, s, λ, τ₀, β, log2n, Δτ, k), std::ios::binary);
  std::ifstream rtfile(logFourierFile("Rt", p, s, λ, τ₀, β, log2n, Δτ, k), std::ios::binary);

  if ((!cfile.is_open() || !rfile.is_open()) || (!ctfile.is_open() || !rtfile.is_open())) {
    return false;
  }

  unsigned N = pow(2, log2n);

  cfile.read((char*)(C.data()), N * sizeof(Real));
  cfile.close();

  rfile.read((char*)(R.data()), N * sizeof(Real));
  rfile.close();

  ctfile.read((char*)(Ct.data()), N * sizeof(Complex));
  ctfile.close();

  rtfile.read((char*)(Rt.data()), N * sizeof(Complex));
  rtfile.close();

  return true;
}

std::tuple<std::vector<Complex>, std::vector<Complex>> RddfCtdfCt(LogarithmicFourierTransform& fft, const std::vector<Real>& C, const std::vector<Real>& R, unsigned p, unsigned s, Real λ) {
  std::vector<Real> dfC(C.size());
  std::vector<Real> RddfC(C.size());
  for (unsigned n = 0; n < C.size(); n++) {
    RddfC[n] = R[n] * ddf(λ, p, s, C[n]);
    dfC[n] = df(λ, p, s, C[n]);
  }
  std::vector<Complex> RddfCt = fft.fourier(RddfC, false);
  std::vector<Complex> dfCt = fft.fourier(dfC, true);

  return {RddfCt, dfCt};
}

Real estimateZ(LogarithmicFourierTransform& fft, const std::vector<Real>& C, const std::vector<Complex>& Ct, const std::vector<Real>& R, const std::vector<Complex>& Rt, unsigned p, unsigned s, Real λ, Real τ₀, Real β) {
  auto [RddfCt, dfCt] = RddfCtdfCt(fft, C, R, p, s, λ);
  Real Γ₀ = 1.0;

  return ((2 * Γ₀ * std::conj(Rt[0]) + pow(β, 2) * (RddfCt[0] * Ct[0] + dfCt[0] * std::conj(Rt[0]))) / Ct[0]).real();
}

Real energy(const LogarithmicFourierTransform& fft, std::vector<Real>&  C, const std::vector<Real>& R, unsigned p, unsigned s, Real λ, Real β) {
  Real E = 0;
  for (unsigned n = 0; n < C.size()/2-1; n++) {
    Real h₂ₙ   = fft.t(2*n+1) - fft.t(2*n);
    Real h₂ₙ₊₁ = fft.t(2*n+2) - fft.t(2*n+1);
    Real f₂ₙ   = R[2*n]   * df(λ, p, s, C[2*n]);
    Real f₂ₙ₊₁ = R[2*n+1] * df(λ, p, s, C[2*n+1]);
    Real f₂ₙ₊₂ = R[2*n+2] * df(λ, p, s, C[2*n+2]);
    E += (h₂ₙ + h₂ₙ₊₁) / 6 * (
          (2 - h₂ₙ₊₁ / h₂ₙ) * f₂ₙ
          + pow(h₂ₙ + h₂ₙ₊₁, 2) / (h₂ₙ * h₂ₙ₊₁) * f₂ₙ₊₁
          + (2 - h₂ₙ / h₂ₙ₊₁) * f₂ₙ₊₂
        );
  }
  return  β * E;
}