diff options
author | Jaron Kent-Dobias <jaron@kent-dobias.com> | 2025-04-03 15:45:29 -0300 |
---|---|---|
committer | Jaron Kent-Dobias <jaron@kent-dobias.com> | 2025-04-03 15:45:29 -0300 |
commit | 8364b9a96b8f60a3effbe8843c89717a89a8fc5e (patch) | |
tree | daf2c354746cf85499ceca67f324325574e59f6a | |
parent | 60e0f9b63c9825a265a93a278bf018d301063ea2 (diff) | |
download | code-8364b9a96b8f60a3effbe8843c89717a89a8fc5e.tar.gz code-8364b9a96b8f60a3effbe8843c89717a89a8fc5e.tar.bz2 code-8364b9a96b8f60a3effbe8843c89717a89a8fc5e.zip |
Split functions into library to deduplicate code
-rw-r--r-- | Makefile | 13 | ||||
-rw-r--r-- | fourier.cpp | 86 | ||||
-rw-r--r-- | fourier.hpp | 32 | ||||
-rw-r--r-- | fourier_integrator.cpp | 121 | ||||
-rw-r--r-- | get_energy.cpp | 101 |
5 files changed, 161 insertions, 192 deletions
@@ -1,6 +1,6 @@ all: walk correlations integrator fourier_integrator get_energy -CC := clang++ -std=c++17 -Wno-mathematical-notation-identifier-extension -O3 -march=native -mtune=native -fopenmp -lfftw3 -g +CC := clang++ -std=c++17 -flto -Wno-mathematical-notation-identifier-extension -O3 -march=native -mtune=native -fopenmp walk: walk.cpp $(CC) walk.cpp -o walk @@ -11,8 +11,11 @@ correlations: correlations.cpp integrator: integrator.cpp $(CC) integrator.cpp -o integrator -fourier_integrator: fourier_integrator.cpp - $(CC) fourier_integrator.cpp -o fourier_integrator +fourier.o: fourier.cpp fourier.hpp + $(CC) fourier.cpp -c -o fourier.o -get_energy: get_energy.cpp - $(CC) get_energy.cpp -o get_energy +fourier_integrator: fourier_integrator.cpp fourier.hpp fourier.o + $(CC) fourier_integrator.cpp fourier.o -lfftw3 -o fourier_integrator + +get_energy: get_energy.cpp fourier.hpp fourier.o + $(CC) get_energy.cpp fourier.o -lfftw3 -o get_energy diff --git a/fourier.cpp b/fourier.cpp new file mode 100644 index 0000000..7808989 --- /dev/null +++ b/fourier.cpp @@ -0,0 +1,86 @@ +#include "fourier.hpp" + +inline Real fP(unsigned p, Real q) { + return 0.5 * pow(q, p); +} + +inline Real dfP(unsigned p, Real q) { + return 0.5 * p * pow(q, p - 1); +} + +inline Real ddfP(unsigned p, Real q) { + return 0.5 * p * (p - 1) * pow(q, p - 2); +} + +Real f(Real λ, unsigned p, unsigned s, Real q) { + return (1 - λ) * fP(p, q) + λ * fP(s, q); +} + +Real df(Real λ, unsigned p, unsigned s, Real q) { + return (1 - λ) * dfP(p, q) + λ * dfP(s, q); +} + +Real ddf(Real λ, unsigned p, unsigned s, Real q) { + return (1 - λ) * ddfP(p, q) + λ * ddfP(s, q); +} + +FourierTransform::FourierTransform(unsigned n, Real Δω, Real Δτ, unsigned flags) : a(2 * n), â(n + 1), Δω(Δω), Δτ(Δτ) { + plan_r2c = fftw_plan_dft_r2c_1d(2 * n, a.data(), reinterpret_cast<fftw_complex*>(â.data()), flags); + plan_c2r = fftw_plan_dft_c2r_1d(2 * n, reinterpret_cast<fftw_complex*>(â.data()), a.data(), flags); +} + +FourierTransform::~FourierTransform() { + fftw_destroy_plan(plan_r2c); + fftw_destroy_plan(plan_c2r); + fftw_cleanup(); +} + +std::vector<Complex> FourierTransform::fourier(const std::vector<Real>& c) { + a = c; + fftw_execute(plan_r2c); + std::vector<Complex> ĉ(â.size()); + for (unsigned i = 0; i < â.size(); i++) { + ĉ[i] = â[i] * (Δτ * M_PI); + } + return ĉ; +} + +std::vector<Real> FourierTransform::inverse(const std::vector<Complex>& ĉ) { + â = ĉ; + fftw_execute(plan_c2r); + std::vector<Real> c(a.size()); + for (unsigned i = 0; i < a.size(); i++) { + c[i] = a[i] * (Δω / (2 * M_PI)); + } + return c; +} + +std::string fourierFile(std::string prefix, unsigned p, unsigned s, Real λ, Real τ₀, Real y, unsigned log2n, Real τₘₐₓ) { + return prefix + "_" + std::to_string(p) + "_" + std::to_string(s) + "_" + std::to_string(λ) + "_" + std::to_string(τ₀) + "_" + std::to_string(y) + "_" + std::to_string(log2n) + "_" + std::to_string(τₘₐₓ) + ".dat"; +} + +Real energy(const std::vector<Real>& C, const std::vector<Real>& R, unsigned p, unsigned s, Real λ, Real y, Real Δτ) { + Real e = 0; + + for (unsigned i = 0; i < C.size() / 2; i++) { + e += y * R[i] * df(λ, p, s, C[i]) * M_PI * Δτ; + } + + return e; +} + +std::tuple<std::vector<Complex>, std::vector<Complex>> RddfCtdfCt(FourierTransform& fft, const std::vector<Real>& C, const std::vector<Real>& R, unsigned p, unsigned s, Real λ) { + std::vector<Real> RddfC(C.size()); + for (unsigned i = 0; i < C.size(); i++) { + RddfC[i] = R[i] * ddf(λ, p, s, C[i]); + } + std::vector<Complex> RddfCt = fft.fourier(RddfC); + + std::vector<Real> dfC(C.size()); + for (unsigned i = 0; i < C.size(); i++) { + dfC[i] = df(λ, p, s, C[i]); + } + std::vector<Complex> dfCt = fft.fourier(dfC); + + return {RddfCt, dfCt}; +} diff --git a/fourier.hpp b/fourier.hpp new file mode 100644 index 0000000..86335f7 --- /dev/null +++ b/fourier.hpp @@ -0,0 +1,32 @@ +#pragma once +#include <cmath> +#include <fftw3.h> +#include <complex> +#include <vector> + +using Real = double; +using Complex = std::complex<Real>; +using namespace std::complex_literals; + +Real f(Real λ, unsigned p, unsigned s, Real q); +Real df(Real λ, unsigned p, unsigned s, Real q); +Real ddf(Real λ, unsigned p, unsigned s, Real q); + +class FourierTransform { +private: + std::vector<Real> a; + std::vector<Complex> â; + fftw_plan plan_r2c; + fftw_plan plan_c2r; + Real Δω; + Real Δτ; +public: + FourierTransform(unsigned n, Real Δω, Real Δτ, unsigned flags = 0); + ~FourierTransform(); + std::vector<Complex> fourier(const std::vector<Real>& c); + std::vector<Real> inverse(const std::vector<Complex>& ĉ); +}; + +std::string fourierFile(std::string prefix, unsigned p, unsigned s, Real λ, Real τ₀, Real y, unsigned log2n, Real τₘₐₓ); +Real energy(const std::vector<Real>& C, const std::vector<Real>& R, unsigned p, unsigned s, Real λ, Real y, Real Δτ); +std::tuple<std::vector<Complex>, std::vector<Complex>> RddfCtdfCt(FourierTransform& fft, const std::vector<Real>& C, const std::vector<Real>& R, unsigned p, unsigned s, Real λ); diff --git a/fourier_integrator.cpp b/fourier_integrator.cpp index 22866f2..86c3b24 100644 --- a/fourier_integrator.cpp +++ b/fourier_integrator.cpp @@ -1,80 +1,8 @@ +#include "fourier.hpp" #include <getopt.h> -#include <vector> -#include <cmath> #include <iostream> -#include <fftw3.h> -#include <complex> #include <fstream> -using Real = double; -using Complex = std::complex<Real>; -using namespace std::complex_literals; - -inline Real fP(unsigned p, Real q) { - return 0.5 * pow(q, p); -} - -inline Real dfP(unsigned p, Real q) { - return 0.5 * p * pow(q, p - 1); -} - -inline Real ddfP(unsigned p, Real q) { - return 0.5 * p * (p - 1) * pow(q, p - 2); -} - -inline Real f(Real λ, unsigned p, unsigned s, Real q) { - return (1 - λ) * fP(p, q) + λ * fP(s, q); -} - -inline Real df(Real λ, unsigned p, unsigned s, Real q) { - return (1 - λ) * dfP(p, q) + λ * dfP(s, q); -} - -inline Real ddf(Real λ, unsigned p, unsigned s, Real q) { - return (1 - λ) * ddfP(p, q) + λ * ddfP(s, q); -} - -class FourierTransform { -private: - std::vector<Real> a; - std::vector<Complex> â; - fftw_plan plan_r2c; - fftw_plan plan_c2r; - Real Δω; - Real Δτ; -public: - FourierTransform(unsigned n, Real Δω, Real Δτ) : a(2 * n), â(n + 1), Δω(Δω), Δτ(Δτ) { - plan_r2c = fftw_plan_dft_r2c_1d(2 * n, a.data(), reinterpret_cast<fftw_complex*>(â.data()), 0); - plan_c2r = fftw_plan_dft_c2r_1d(2 * n, reinterpret_cast<fftw_complex*>(â.data()), a.data(), 0); - } - - ~FourierTransform() { - fftw_destroy_plan(plan_r2c); - fftw_destroy_plan(plan_c2r); - fftw_cleanup(); - } - - std::vector<Complex> fourier(const std::vector<Real>& c) { - a = c; - fftw_execute(plan_r2c); - std::vector<Complex> ĉ(â.size()); - for (unsigned i = 0; i < â.size(); i++) { - ĉ[i] = â[i] * (Δτ * M_PI); - } - return ĉ; - } - - std::vector<Real> inverse(const std::vector<Complex>& ĉ) { - â = ĉ; - fftw_execute(plan_c2r); - std::vector<Real> c(a.size()); - for (unsigned i = 0; i < a.size(); i++) { - c[i] = a[i] * (Δω / (2 * M_PI)); - } - return c; - } -}; - int main(int argc, char* argv[]) { unsigned p = 2; unsigned s = 2; @@ -147,6 +75,10 @@ int main(int argc, char* argv[]) { std::vector<Real> R(2 * n); FourierTransform fft(n, Δω, Δτ); + std::vector<Complex> Ct; + std::vector<Complex> Rt; + + Real y = y₀; if (!loadData) { // start from the exact solution for τ₀ = 0 @@ -158,37 +90,30 @@ int main(int argc, char* argv[]) { } R[i] = exp(-z * τ); } + Ct = fft.fourier(C); + Rt = fft.fourier(R); } else { - std::string file_end = std::to_string(p) + "_" + std::to_string(s) + "_" + std::to_string(λ) + "_" + std::to_string(τ₀) + "_" + std::to_string(y₀) + "_" + std::to_string(log2n) + "_" + std::to_string(τₘₐₓ) + ".dat"; - std::ifstream cfile("C_"+file_end, std::ios::binary); + std::ifstream cfile(fourierFile("C", p, s, λ, τ₀, y, log2n, τₘₐₓ), std::ios::binary); cfile.read((char*)(C.data()), C.size() * sizeof(Real)); cfile.close(); - std::ifstream rfile("R_"+file_end, std::ios::binary); + std::ifstream rfile(fourierFile("R", p, s, λ, τ₀, y, log2n, τₘₐₓ), std::ios::binary); rfile.read((char*)(R.data()), R.size() * sizeof(Real)); rfile.close(); - } - std::vector<Complex> Ct = fft.fourier(C); - std::vector<Complex> Rt = fft.fourier(R); + Ct = fft.fourier(C); + Rt = fft.fourier(R); - Real y = y₀; + auto [RddfCt, dfCt] = RddfCtdfCt(fft, C, R, p, s, λ); + + z = ((Γ₀ * std::conj(Rt[0]) + pow(y, 2) * (RddfCt[0] * Ct[0] + dfCt[0] * std::conj(Rt[0]))) / Ct[0]).real(); + } while (y += Δy, y <= yₘₐₓ) { Real ΔC = 1;; unsigned it = 0; while (sqrt(ΔC / C.size()) > ε) { it++; - std::vector<Real> RddfC(C.size()); - for (unsigned i = 0; i < C.size(); i++) { - RddfC[i] = R[i] * ddf(λ, p, s, C[i]); - } - std::vector<Complex> RddfCt = fft.fourier(RddfC); - - std::vector<Real> dfC(C.size()); - for (unsigned i = 0; i < C.size(); i++) { - dfC[i] = df(λ, p, s, C[i]); - } - std::vector<Complex> dfCt = fft.fourier(dfC); + auto [RddfCt, dfCt] = RddfCtdfCt(fft, C, R, p, s, λ); for (unsigned i = 0; i < Rt.size(); i++) { Real ω = i * Δω; @@ -230,21 +155,15 @@ int main(int argc, char* argv[]) { } } - Real energy = 0; - - for (unsigned i = 0; i < n; i++) { - energy += y * R[i] * df(λ, p, s, C[i]) * M_PI * Δτ; - } - - std::cerr << "y " << y << " " << energy << std::endl; + Real e = energy(C, R, p, s, λ, y, Δτ); - std::string file_end = std::to_string(p) + "_" + std::to_string(s) + "_" + std::to_string(λ) + "_" + std::to_string(τ₀) + "_" + std::to_string(y) + "_" + std::to_string(log2n) + "_" + std::to_string(τₘₐₓ) + ".dat"; + std::cerr << "y " << y << " " << e << " " << z << std::endl; - std::ofstream outfile("C_" + file_end, std::ios::out | std::ios::binary); + std::ofstream outfile(fourierFile("C", p, s, λ, τ₀, y, log2n, τₘₐₓ), std::ios::out | std::ios::binary); outfile.write((const char*)(C.data()), C.size() * sizeof(Real)); outfile.close(); - std::ofstream outfileR("R_" + file_end, std::ios::out | std::ios::binary); + std::ofstream outfileR(fourierFile("R", p, s, λ, τ₀, y, log2n, τₘₐₓ), std::ios::out | std::ios::binary); outfileR.write((const char*)(R.data()), R.size() * sizeof(Real)); outfileR.close(); } diff --git a/get_energy.cpp b/get_energy.cpp index 6164a6d..a388e0c 100644 --- a/get_energy.cpp +++ b/get_energy.cpp @@ -1,80 +1,8 @@ +#include "fourier.hpp" #include <getopt.h> -#include <vector> -#include <cmath> #include <iostream> -#include <fftw3.h> -#include <complex> #include <fstream> -using Real = double; -using Complex = std::complex<Real>; -using namespace std::complex_literals; - -inline Real fP(unsigned p, Real q) { - return 0.5 * pow(q, p); -} - -inline Real dfP(unsigned p, Real q) { - return 0.5 * p * pow(q, p - 1); -} - -inline Real ddfP(unsigned p, Real q) { - return 0.5 * p * (p - 1) * pow(q, p - 2); -} - -inline Real f(Real λ, unsigned p, unsigned s, Real q) { - return (1 - λ) * fP(p, q) + λ * fP(s, q); -} - -inline Real df(Real λ, unsigned p, unsigned s, Real q) { - return (1 - λ) * dfP(p, q) + λ * dfP(s, q); -} - -inline Real ddf(Real λ, unsigned p, unsigned s, Real q) { - return (1 - λ) * ddfP(p, q) + λ * ddfP(s, q); -} - -class FourierTransform { -private: - std::vector<Real> a; - std::vector<Complex> â; - fftw_plan plan_r2c; - fftw_plan plan_c2r; - Real Δω; - Real Δτ; -public: - FourierTransform(unsigned n, Real Δω, Real Δτ, unsigned flags = 0) : a(2 * n), â(n + 1), Δω(Δω), Δτ(Δτ) { - plan_r2c = fftw_plan_dft_r2c_1d(2 * n, a.data(), reinterpret_cast<fftw_complex*>(â.data()), flags); - plan_c2r = fftw_plan_dft_c2r_1d(2 * n, reinterpret_cast<fftw_complex*>(â.data()), a.data(), flags); - } - - ~FourierTransform() { - fftw_destroy_plan(plan_r2c); - fftw_destroy_plan(plan_c2r); - fftw_cleanup(); - } - - std::vector<Complex> fourier(const std::vector<Real>& c) { - a = c; - fftw_execute(plan_r2c); - std::vector<Complex> ĉ(â.size()); - for (unsigned i = 0; i < â.size(); i++) { - ĉ[i] = â[i] * (Δτ * M_PI); - } - return ĉ; - } - - std::vector<Real> inverse(const std::vector<Complex>& ĉ) { - â = ĉ; - fftw_execute(plan_c2r); - std::vector<Real> c(a.size()); - for (unsigned i = 0; i < a.size(); i++) { - c[i] = a[i] * (Δω / (2 * M_PI)); - } - return c; - } -}; - int main(int argc, char* argv[]) { unsigned p = 2; unsigned s = 2; @@ -133,24 +61,25 @@ int main(int argc, char* argv[]) { std::vector<Real> C(2 * n); std::vector<Real> R(2 * n); - std::string file_end = std::to_string(p) + "_" + std::to_string(s) + "_" + std::to_string(λ) + "_" + std::to_string(τ₀) + "_" + std::to_string(y) + "_" + std::to_string(log2n) + "_" + std::to_string(τₘₐₓ) + ".dat"; - std::ifstream cfile("C_"+file_end, std::ios::binary); + std::ifstream cfile(fourierFile("C", p, s, λ, τ₀, y, log2n, τₘₐₓ), std::ios::binary); if (cfile.is_open()) { - cfile.read((char*)(C.data()), C.size() * sizeof(Real)); - cfile.close(); - std::ifstream rfile("R_"+file_end, std::ios::binary); - rfile.read((char*)(R.data()), R.size() * sizeof(Real)); - rfile.close(); + cfile.read((char*)(C.data()), C.size() * sizeof(Real)); + cfile.close(); - Real energy = 0; + std::ifstream rfile(fourierFile("R", p, s, λ, τ₀, y, log2n, τₘₐₓ), std::ios::binary); + rfile.read((char*)(R.data()), R.size() * sizeof(Real)); + rfile.close(); - for (unsigned i = 0; i < n; i++) { - energy += y * R[i] * df(λ, p, s, C[i]) * M_PI * Δτ; - } + Real e = energy(C, R, p, s, λ, y, Δτ); + + std::vector<Complex> Ct = fft.fourier(C); + std::vector<Complex> Rt = fft.fourier(R); + + auto [RddfCt, dfCt] = RddfCtdfCt(fft, C, R, p, s, λ); - std::vector<Complex> Ct = fft.fourier(C); + Real z = ((std::conj(Rt[0]) + pow(y, 2) * (RddfCt[0] * Ct[0] + dfCt[0] * std::conj(Rt[0]))) / Ct[0]).real(); - std::cout << y << " " << energy << " " << Ct[0].real() << std::endl; + std::cout << y << " " << e << " " << Ct[0].real() << " " << z << std::endl; } } |