Triumvirate C++ API 0.5.0
Three-point clustering measurements in large-scale structure analyses.
Loading...
Searching...
No Matches
maths.cpp
Go to the documentation of this file.
1// Copyright (C) [GPLv3 Licence]
2//
3// This file is part of the Triumvirate program. See the COPYRIGHT
4// and LICENCE files at the top-level directory of this distribution
5// for details of copyright and licensing.
6//
7// This program is free software: you can redistribute it and/or modify it
8// under the terms of the GNU General Public License as published by
9// the Free Software Foundation, either version 3 of the License, or
10// (at your option) any later version.
11//
12// This program is distributed in the hope that it will be useful,
13// but WITHOUT ANY WARRANTY; without even the implied warranty of
14// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
15// See the GNU General Public License for more details.
16//
17// You should have received a copy of the GNU General Public License
18// along with this program. If not, see <https://www.gnu.org/licenses/>.
19
27#include "maths.hpp"
28
29namespace trvs = trv::sys;
30
31namespace trv {
32namespace maths {
33
34// ***********************************************************************
35// Complex numbers
36// ***********************************************************************
37
39const std::complex<double> M_I(0., 1.);
41
42std::complex<double> eval_complex_in_polar(double r, double theta) {
43 return r * (std::cos(theta) + M_I * std::sin(theta));
44}
45
46
47// ***********************************************************************
48// Vectors
49// ***********************************************************************
50
51double get_vec3d_magnitude(std::vector<double> vec) {
52 return std::sqrt(vec[0] * vec[0] + vec[1] * vec[1] + vec[2] * vec[2]);
53}
54
55double get_vec3d_magnitude(double* vec) {
56 return std::sqrt(vec[0] * vec[0] + vec[1] * vec[1] + vec[2] * vec[2]);
57}
58
59
60// ***********************************************************************
61// Gamma function
62// ***********************************************************************
63
64// -----------------------------------------------------------------------
65// Lanzcos approximation
66// -----------------------------------------------------------------------
67
69const double gconst_lanczos = 7.;
70
71std::complex<double> eval_lanczos_approx_series(std::complex<double> z) {
72 // CAVEAT: Discretionary choices (for the Lanczos approximation series).
73 // Declare the number of approximation terms and the corresponding
74 // Lanczos approximation coefficients.
75 const int nterm_lanczos = 9;
76 const double pcoeff_lanczos[] = {
77 0.999999999999809932277,
78 676.520368121885098567,
79 -1259.13921672240287047,
80 771.323428777653078849,
81 -176.615029162140599066,
82 12.5073432786869048145,
83 -0.138571095265720116896,
84 9.98436957801957085956e-6,
85 1.50563273514931155834e-7,
86 }; // NOTE: constant number of significant figures.
87
88 std::complex<double> series = pcoeff_lanczos[0];
89 for (int i = 1; i < nterm_lanczos; i++) {
90 series += pcoeff_lanczos[i] / (z + double(i));
91 }
92
93 return series;
94}
95
96std::complex<double> eval_gamma_lanczos(std::complex<double> z) {
97 // Exploit Euler's reflection formula as the Lanczos approximation
98 // is only valid for Re{z} > 1/2.
99 if (z.real() < 1./2) {
100 return M_PI / (std::sin(M_PI * z) * eval_gamma_lanczos(1. - z));
101 }
102
103 // Substitute variables into the approximation formula.
104 z -= 1.;
105
106 std::complex<double> t = z + gconst_lanczos + 1./2;
107 std::complex<double> series = eval_lanczos_approx_series(z);
108
109 std::complex<double> gamma = std::sqrt(2.*M_PI)
110 * std::pow(t, z + 1./2) * std::exp(-t) * series;
111
112 return gamma;
113}
114
115
116// -----------------------------------------------------------------------
117// Component evaluation
118// -----------------------------------------------------------------------
119
120void get_lngamma_parts(double x, double y, double& lnr, double& theta) {
121 gsl_sf_result lnr_result, theta_result;
122 gsl_sf_lngamma_complex_e(x, y, &lnr_result, &theta_result);
123
124 lnr = lnr_result.val;
125 theta = theta_result.val;
126}
127
128std::complex<double> eval_gamma_lnratio(double mu, std::complex<double> nu) {
129 std::complex<double> x_p = (mu + 1. + nu)/2.;
130 std::complex<double> x_m = (mu + 1. - nu)/2.;
131
132 // Although the Stirling approximation is used, the Lanczos
133 // approximation for log-gamma should remain accurate enough.
134 const double cutoff_asymp = 100.;
135 std::complex<double> lnratio;
136 if (nu.imag() < cutoff_asymp) {
137 double lnr_p, theta_p;
138 double lnr_m, theta_m;
139
140 get_lngamma_parts(x_p.real(), x_p.imag(), lnr_p, theta_p);
141 get_lngamma_parts(x_m.real(), x_m.imag(), lnr_m, theta_m);
142
143 lnratio.real(lnr_p - lnr_m);
144 lnratio.imag(theta_p - theta_m);
145 } else {
146 lnratio =
147 - nu
148 + (x_p - 1./2) * std::log(x_p) - (x_m - 1./2) * std::log(x_m)
149 + 1./12 * (1./x_p - 1./x_m)
150 - 1./360 * (1./std::pow(x_p, 3) - 1./std::pow(x_m, 3))
151 + 1./1260 * (1./std::pow(x_p, 5) - 1./std::pow(x_m, 5));
152 }
153
154 return lnratio;
155}
156
157
158// ***********************************************************************
159// Spherical harmonics
160// ***********************************************************************
161
162// CAVEAT: Discretionary choice such that eps = 1.e-9.
163const double eps_coupling = 1.e-9;
164
165double wigner_3j(int j1, int j2, int j3, int m1, int m2, int m3) {
166 return gsl_sf_coupling_3j(2*j1, 2*j2, 2*j3, 2*m1, 2*m2, 2*m3);
167}
168
169std::complex<double> \
170SphericalHarmonicCalculator::calc_reduced_spherical_harmonic(
171 const int ell, const int m, double pos[3]
172) {
173 // CAVEAT: Discretionary choice such that eps = 1.e-9.
174 const double eps = 1.e-9;
175
176 // Return unity in the trivial case.
177 if (ell == 0 && m == 0) {return 1.;}
178
179 // Calculate modulus.
180 double xyz_mod_sq = 0.;
181 for (int iaxis = 0; iaxis < 3; iaxis++) {
182 xyz_mod_sq += pos[iaxis] * pos[iaxis];
183 } // r² = x² + y² + z²
184
185 double xyz_mod = std::sqrt(xyz_mod_sq); // r = √(x² + y² + z²)
186
187 // Return zero in the trivial case.
188 if (std::fabs(xyz_mod) < eps) {return 0.;}
189
190 // Calculate the angular variable μ = cos(θ).
191 double mu = pos[2] / xyz_mod; // μ = z / r
192
193 // Calculate the angular variable Ï•.
194 double xy_mod = std::sqrt(pos[0] * pos[0] + pos[1] * pos[1]); // r_xy =
195 // √(x² + y²)
196
197 double phi = 0.;
198 if (std::fabs(xy_mod) >= eps) {
199 phi = std::acos(pos[0] / xy_mod); // Ï• = arccos(x / r_xy)
200 if (pos[1] < 0.) {
201 phi = - phi + 2.*M_PI; // ϕ ∈ [π, 2π] if y < 0
202 }
203 }
204
205 // Calculate spherical harmonics with m >= 0 via the normalised
206 // associated Legendre polynomial, i.e. Y_lm = √((2l + 1)/(4π))
207 // √((l - |m|)!/(l + |m|)!) P_l^|m|(μ) * exp(imϕ).
208 std::complex<double> ylm = std::exp(M_I * double(m) * phi)
209 * gsl_sf_legendre_sphPlm(ell, std::abs(m), mu);
210
211 // Impose parity and conjugation.
212 ylm = std::pow(-1, (m - std::abs(m))/2) * std::conj(ylm);
213
214 // Normalise to the reduced form.
215 ylm *= std::sqrt(4.*M_PI / (2.*ell + 1.));
216
217 return ylm;
218}
219
220void SphericalHarmonicCalculator::\
221store_reduced_spherical_harmonic_in_fourier_space(
222 const int ell, const int m,
223 const double boxsize[3], const int ngrid[3],
224 std::vector< std::complex<double> >& ylm_out
225) {
226 // Determine the fundamental wavenumber in each dimension.
227 double dk[3] = {
228 2.*M_PI / boxsize[0], 2.*M_PI / boxsize[1], 2.*M_PI / boxsize[2]
229 };
230
231 // Assign a wavevector to each grid cell.
232#ifdef TRV_USE_OMP
233#pragma omp parallel for collapse(3)
234#endif // TRV_USE_OMP
235 for (int i = 0; i < ngrid[0]; i++) {
236 for (int j = 0; j < ngrid[1]; j++) {
237 for (int k = 0; k < ngrid[2]; k++) {
238 // Lay the 'bricks' vertically, then inwards, then to
239 // the right, i.e. along z-axis, y-axis and then x-axis.
240 // The assigned flattened-grid array index is
241 // (i * ngrid_y * ngrid_z + j * ngrid_z + k)
242 // where ngrid is the grid number along each axis.
243 long long idx_grid =
244 (i * static_cast<long long>(ngrid[1]) + j) * ngrid[2] + k;
245
246 // This conforms to the (absurd) FFT array-ordering convention
247 // that negative wavenumbers/frequencies come after zero and
248 // positive wavenumbers/frequencies.
249 double kvec[3];
250 kvec[0] = (i < ngrid[0]/2) ? i * dk[0] : (i - ngrid[0]) * dk[0];
251 kvec[1] = (j < ngrid[1]/2) ? j * dk[1] : (j - ngrid[1]) * dk[1];
252 kvec[2] = (k < ngrid[2]/2) ? k * dk[2] : (k - ngrid[2]) * dk[2];
253
254 ylm_out[idx_grid] = calc_reduced_spherical_harmonic(ell, m, kvec);
255 }
256 }
257 }
258}
259
260void SphericalHarmonicCalculator::\
261store_reduced_spherical_harmonic_in_config_space(
262 const int ell, const int m,
263 const double boxsize[3], const int ngrid[3],
264 std::vector< std::complex<double> >& ylm_out
265) {
266 // Determine the grid cell size in each dimension.
267 double dr[3] = {
268 boxsize[0] / double(ngrid[0]),
269 boxsize[1] / double(ngrid[1]),
270 boxsize[2] / double(ngrid[2])
271 };
272
273 // Assign a position vector to each grid cell.
274#ifdef TRV_USE_OMP
275#pragma omp parallel for collapse(3)
276#endif // TRV_USE_OMP
277 for (int i = 0; i < ngrid[0]; i++) {
278 for (int j = 0; j < ngrid[1]; j++) {
279 for (int k = 0; k < ngrid[2]; k++) {
280 // Lay the 'bricks' vertically, then inwards, then to
281 // the right, i.e. along z-axis, y-axis and then x-axis.
282 // The assigned flattened-grid array index is
283 // (i * ngrid_y * ngrid_z + j * ngrid_z + k)
284 // where ngrid is the grid number along each axis.
285 long long idx_grid =
286 (i * static_cast<long long>(ngrid[1]) + j) * ngrid[2] + k;
287
288 // This conforms to the (absurd) FFT array-ordering convention
289 // that negative wavenumbers/frequencies come after zero and
290 // positive wavenumbers/frequencies.
291 double rvec[3];
292 rvec[0] = (i < ngrid[0]/2) ? i * dr[0] : (i - ngrid[0]) * dr[0];
293 rvec[1] = (j < ngrid[1]/2) ? j * dr[1] : (j - ngrid[1]) * dr[1];
294 rvec[2] = (k < ngrid[2]/2) ? k * dr[2] : (k - ngrid[2]) * dr[2];
295
296 ylm_out[idx_grid] = calc_reduced_spherical_harmonic(ell, m, rvec);
297 }
298 }
299 }
300}
301
302
303// ***********************************************************************
304// Spherical Bessel function
305// ***********************************************************************
306
308 const int ell
309) : order(ell) {
310 // Set up sampling range and number.
311 this->split = (this->split >= this->order * this->order) ?
312 this->split : this->order * this->order;
313
314 const double xmin = 0.; // minimum of interpolation range
315 const double xmax = this->split; // maximum of interpolation range
316 const double dx = this->step; // interpolation step size
317
318 int nsample = int((xmax - xmin)/dx) + 1; // interpolation sample number
319
320 // Initialise and evaluate at sample points.
321 double* x = new double[nsample];
322 double* j_ell = new double[nsample];
323
324#ifdef TRV_USE_OMP
325#pragma omp parallel for
326#endif // TRV_USE_OMP
327 for (int i = 0; i < nsample; i++) {
328 x[i] = xmin + dx * i;
329 j_ell[i] = gsl_sf_bessel_jl(this->order, x[i]);
330 }
331
332 // Initialise the interpolator using cubic spline and the accelerator.
333 this->accel = gsl_interp_accel_alloc();
334 this->spline = gsl_spline_alloc(gsl_interp_cspline, nsample);
335
336 gsl_spline_init(this->spline, x, j_ell, nsample);
337
338 delete[] x; delete[] j_ell;
339}
340
342 const SphericalBesselCalculator& other
343) {
344 this->order = other.order;
345 this->split = other.split;
346 this->step = other.step;
347
348 this->accel = gsl_interp_accel_alloc();
349 this->spline = gsl_spline_alloc(gsl_interp_cspline, other.spline->size);
350
351 gsl_spline_init(
352 this->spline, other.spline->x, other.spline->y, other.spline->size
353 );
354}
355
357 if (this->accel != nullptr) {
358 gsl_interp_accel_free(this->accel); this->accel = nullptr;
359 }
360
361 if (this->spline != nullptr) {
362 gsl_spline_free(this->spline); this->spline = nullptr;
363 }
364}
365
367 if (x >= this->split) {
368 return gsl_sf_bessel_jl(this->order, x);
369 } else {
370 // NOTE: This is a computational bottleneck.
371 return gsl_spline_eval(this->spline, x, this->accel);
372 }
373}
374
375} // namespace trv::maths
376} // namespace trv
Interpolated spherical Bessel function of the first kind.
Definition maths.hpp:258
double eval(double x)
Evaluate the interpolated function.
Definition maths.cpp:366
SphericalBesselCalculator(const int ell)
Construct the interpolated function calculator.
Definition maths.cpp:307
~SphericalBesselCalculator()
Destruct the interpolated function.
Definition maths.cpp:356
static std::complex< double > calc_reduced_spherical_harmonic(const int ell, const int m, double pos[3])
Calculate the reduced spherical harmonic.
Definition maths.cpp:170
Mathematical calculations.
std::complex< double > eval_complex_in_polar(double r, double theta)
Evaluate a complex number in the polar form.
Definition maths.cpp:42
double get_vec3d_magnitude(std::vector< double > vec)
Return the magnitude of a 3-d vector.
Definition maths.cpp:51
std::complex< double > eval_gamma_lanczos(std::complex< double > z)
Evaluate the gamma function on the complex plane using the Lanczos approximation.
Definition maths.cpp:96
void get_lngamma_parts(double x, double y, double &lnr, double &theta)
Get the real and imaginary parts of the log-gamma function .
Definition maths.cpp:120
std::complex< double > eval_gamma_lnratio(double mu, std::complex< double > nu)
Evaluate the logarithm of the ratio of two gamma functions.
Definition maths.cpp:128
std::complex< double > eval_lanczos_approx_series(std::complex< double > z)
Evaluate the Lanczos approximation series for the gamma function.
Definition maths.cpp:71
const std::complex< double > M_I
imaginary unit
const double eps_coupling
zero-tolerance for Wigner 3-j coupling coefficients
Definition maths.cpp:163
const double gconst_lanczos
Lanczos approximation constant.
Definition maths.cpp:69
double wigner_3j(int j1, int j2, int j3, int m1, int m2, int m3)
Calculate Wigner 3-j symbol.
Definition maths.cpp:165