31const double eps_norm = 1.e-5;
44 return (2*ell + 1) * (2*ELL + 1)
57 if (particles.
pdata ==
nullptr) {
67#if defined(__GNUC__) && !defined(__clang__)
68#pragma omp parallel for simd reduction(+:norm)
70#pragma omp parallel for reduction(+:norm)
73 for (
int pid = 0; pid < particles.
ntotal; pid++) {
74 norm += particles[pid].ws
75 * particles[pid].nz * std::pow(particles[pid].wc, 2);
81 "Particle 'nz' values appear to be all zeros. "
82 "Check the input catalogue contains valid 'nz' field."
86 "Particle 'nz' values appear to be all zeros. "
87 "Check the input catalogue contains valid 'nz' field.\n"
91 double norm_factor = 1. / (alpha * norm);
104 norm_factor /= std::pow(alpha, 2);
118 fftw_complex* weight_data =
nullptr;
119 weight_data = fftw_alloc_complex(particles_data.
ntotal);
121 fftw_complex* weight_rand =
nullptr;
122 weight_rand = fftw_alloc_complex(particles_rand.
ntotal);
129#if defined(__GNUC__) && !defined(__clang__)
130#pragma omp parallel for simd
132#pragma omp parallel for
135 for (
int pid = 0; pid < particles_data.
ntotal; pid++) {
136 weight_data[pid][0] = particles_data[pid].w;
137 weight_data[pid][1] = 0.;
141#if defined(__GNUC__) && !defined(__clang__)
142#pragma omp parallel for simd
144#pragma omp parallel for
147 for (
int pid = 0; pid < particles_rand.
ntotal; pid++) {
148 weight_rand[pid][0] = particles_rand[pid].w;
149 weight_rand[pid][1] = 0.;
155 fftw_free(weight_data); weight_data =
nullptr;
156 fftw_free(weight_rand); weight_rand =
nullptr;
165#if defined(__GNUC__) && !defined(__clang__)
166#pragma omp parallel for simd reduction(+:norm)
168#pragma omp parallel for reduction(+:norm)
171 for (
long long gid = 0; gid < params.
nmesh; gid++) {
172 norm += mesh_data.
field[gid][0] * mesh_rand.
field[gid][0];
175 double vol_cell = params.
volume / double(params.
nmesh);
177 double norm_factor = 1. / (alpha * vol_cell * norm);
186 double padding,
double cellsize,
const std::string& assignment
191 double boxsize_norm = (1. + padding) * std::max(
196 int ngrid_norm = std::ceil(boxsize_norm / cellsize);
198 ngrid_norm += ngrid_norm % 2;
199 boxsize_norm = ngrid_norm * cellsize;
201 for (
int iaxis = 0; iaxis < 3; iaxis++) {
202 params_norm.
boxsize[iaxis] = boxsize_norm;
203 params_norm.
ngrid[iaxis] = ngrid_norm;
212 particles_data, particles_rand, params_norm, alpha
226 if (particles.
pdata ==
nullptr) {
233 double shotnoise = 0.;
236#if defined(__GNUC__) && !defined(__clang__)
237#pragma omp parallel for simd reduction(+:shotnoise)
239#pragma omp parallel for reduction(+:shotnoise)
242 for (
int pid = 0; pid < particles.
ntotal; pid++) {
244 std::pow(particles[pid].ws, 2) * std::pow(particles[pid].wc, 2);
253 double alpha,
int ell,
int m
255 double sn_data_real = 0., sn_data_imag = 0.;
258#pragma omp parallel for reduction(+:sn_data_real, sn_data_imag)
260 for (
int pid = 0; pid < particles_data.
ntotal; pid++) {
262 los_data[pid].
pos[0], los_data[pid].
pos[1], los_data[pid].
pos[2]
265 std::complex<double> ylm = trvm::SphericalHarmonicCalculator::
266 calc_reduced_spherical_harmonic(ell, m, los_);
268 std::complex<double> sn_part = ylm * std::pow(particles_data[pid].w, 2);
269 double sn_part_real = sn_part.real();
270 double sn_part_imag = sn_part.imag();
272 sn_data_real += sn_part_real;
273 sn_data_imag += sn_part_imag;
276 std::complex<double> sn_data(sn_data_real, sn_data_imag);
278 double sn_rand_real = 0., sn_rand_imag = 0.;
281#pragma omp parallel for reduction(+:sn_rand_real, sn_rand_imag)
283 for (
int pid = 0; pid < particles_rand.
ntotal; pid++) {
285 los_rand[pid].
pos[0], los_rand[pid].
pos[1], los_rand[pid].
pos[2]
288 std::complex<double> ylm = trvm::SphericalHarmonicCalculator::
289 calc_reduced_spherical_harmonic(ell, m, los_);
291 std::complex<double> sn_part = ylm * std::pow(particles_rand[pid].w, 2);
292 double sn_part_real = sn_part.real();
293 double sn_part_imag = sn_part.imag();
295 sn_rand_real += sn_part_real;
296 sn_rand_imag += sn_part_imag;
299 std::complex<double> sn_rand(sn_rand_real, sn_rand_imag);
301 return sn_data + std::pow(alpha, 2) * sn_rand;
306 double alpha,
int ell,
int m
308 double sn_real = 0., sn_imag = 0.;
311#pragma omp parallel for reduction(+:sn_real, sn_imag)
313 for (
int pid = 0; pid < particles.
ntotal; pid++) {
314 double los_[3] = {los[pid].
pos[0], los[pid].
pos[1], los[pid].
pos[2]};
316 std::complex<double> ylm = trvm::SphericalHarmonicCalculator::
317 calc_reduced_spherical_harmonic(ell, m, los_);
319 std::complex<double> sn_part = ylm * std::pow(particles[pid].w, 2);
320 double sn_part_real = sn_part.real();
321 double sn_part_imag = sn_part.imag();
323 sn_real += sn_part_real;
324 sn_imag += sn_part_imag;
327 std::complex<double> sn(sn_real, sn_imag);
329 return std::pow(alpha, 2) * sn;
350 "Computing power spectrum from paired survey-type catalogues..."
360 int ell1 = params.
ELL;
363 int* nmodes_save =
new int[kbinning.
num_bins];
364 double* k_save =
new double[kbinning.
num_bins];
365 std::complex<double>* pk_save =
new std::complex<double>[kbinning.
num_bins];
366 std::complex<double>* sn_save =
new std::complex<double>[kbinning.
num_bins];
367 for (
int ibin = 0; ibin < kbinning.
num_bins; ibin++) {
368 nmodes_save[ibin] = 0;
378#if defined(TRV_USE_OMP) && defined(TRV_USE_FFTWOMP)
382 MeshField dn_00(params,
true,
"`dn_00`");
384 catalogue_data, catalogue_rand, los_data, los_rand, alpha, 0, 0
390 for (
int M_ = - params.
ELL; M_ <= params.
ELL; M_++) {
391 MeshField dn_LM(params,
true,
"`dn_LM`");
393 catalogue_data, catalogue_rand, los_data, los_rand, alpha, params.
ELL, M_
398 catalogue_data, catalogue_rand, los_data, los_rand, alpha, params.
ELL, M_
403 for (
int m1 = - ell1; m1 <= ell1; m1++) {
408 dn_LM, dn_00, sn_amp, ell1, m1, kbinning
411 for (
int ibin = 0; ibin < kbinning.
num_bins; ibin++) {
412 pk_save[ibin] += coupling * stats_2pt.
pk[ibin];
413 sn_save[ibin] += coupling * stats_2pt.
sn[ibin];
416 if (M_ == 0 && m1 == 0) {
417 for (
int ibin = 0; ibin < kbinning.
num_bins; ibin++) {
418 nmodes_save[ibin] = stats_2pt.
nmodes[ibin];
419 k_save[ibin] = stats_2pt.
k[ibin];
434 for (
int ibin = 0; ibin < kbinning.
num_bins; ibin++) {
436 powspec_out.
keff.push_back(k_save[ibin]);
437 powspec_out.
nmodes.push_back(nmodes_save[ibin]);
438 powspec_out.
pk_raw.push_back(norm_factor * pk_save[ibin]);
439 powspec_out.
pk_shot.push_back(norm_factor * sn_save[ibin]);
443 delete[] nmodes_save;
delete[] k_save;
delete[] pk_save;
delete[] sn_save;
447 "... computed power spectrum from paired survey-type catalogues."
464 "Computing two-point correlation function from "
465 "paired survey-type catalogues..."
475 int ell1 = params.
ELL;
478 int* npairs_save =
new int[rbinning.
num_bins];
479 double* r_save =
new double[rbinning.
num_bins];
480 std::complex<double>* xi_save =
new std::complex<double>[rbinning.
num_bins];
481 for (
int ibin = 0; ibin < rbinning.
num_bins; ibin++) {
482 npairs_save[ibin] = 0;
491#if defined(TRV_USE_OMP) && defined(TRV_USE_FFTWOMP)
495 MeshField dn_00(params,
true,
"`dn_00`");
497 catalogue_data, catalogue_rand, los_data, los_rand, alpha, 0, 0
503 for (
int M_ = - params.
ELL; M_ <= params.
ELL; M_++) {
504 MeshField dn_LM(params,
true,
"`dn_LM`");
506 catalogue_data, catalogue_rand, los_data, los_rand, alpha, params.
ELL, M_
511 catalogue_data, catalogue_rand, los_data, los_rand, alpha, params.
ELL, M_
516 for (
int m1 = - ell1; m1 <= ell1; m1++) {
521 dn_LM, dn_00, sn_amp, ell1, m1, rbinning
524 for (
int ibin = 0; ibin < rbinning.
num_bins; ibin++) {
525 xi_save[ibin] += coupling * stats_2pt.
xi[ibin];
528 if (M_ == 0 && m1 == 0) {
529 for (
int ibin = 0; ibin < rbinning.
num_bins; ibin++) {
530 npairs_save[ibin] = stats_2pt.
npairs[ibin];
531 r_save[ibin] = stats_2pt.
r[ibin];
538 "Two-point correlation function term computed at order M = %d.", M_
548 for (
int ibin = 0; ibin < rbinning.
num_bins; ibin++) {
550 corrfunc_out.
reff.push_back(r_save[ibin]);
551 corrfunc_out.
npairs.push_back(npairs_save[ibin]);
552 corrfunc_out.
xi.push_back(norm_factor * xi_save[ibin]);
556 delete[] npairs_save;
delete[] r_save;
delete[] xi_save;
560 "... computed two-point correlation function "
561 "from paired survey-type catalogues."
577 "Computing power spectrum from a periodic-box simulation-type catalogue "
578 "in the global plane-parallel approximation."
586 int* nmodes_save =
new int[kbinning.
num_bins];
587 double* k_save =
new double[kbinning.
num_bins];
588 std::complex<double>* pk_save =
new std::complex<double>[kbinning.
num_bins];
589 std::complex<double>* sn_save =
new std::complex<double>[kbinning.
num_bins];
590 for (
int ibin = 0; ibin < kbinning.
num_bins; ibin++) {
591 nmodes_save[ibin] = 0;
598 double norm = double(catalogue_data.
ntotal) * double(catalogue_data.
ntotal)
600 if (std::fabs(1 - norm * norm_factor) > eps_norm) {
603 "Power spectrum normalisation input differs from "
604 "expected value for an unweight field in a periodic box."
613#if defined(TRV_USE_OMP) && defined(TRV_USE_FFTWOMP)
622 std::complex<double> sn_amp = double(catalogue_data.
ntotal);
628 dn, dn, sn_amp, params.
ELL, 0, kbinning
631 for (
int ibin = 0; ibin < kbinning.
num_bins; ibin++) {
632 nmodes_save[ibin] = stats_2pt.
nmodes[ibin];
633 k_save[ibin] = stats_2pt.
k[ibin];
634 pk_save[ibin] += double(2*params.
ELL + 1) * stats_2pt.
pk[ibin];
635 sn_save[ibin] += double(2*params.
ELL + 1) * stats_2pt.
sn[ibin];
644 for (
int ibin = 0; ibin < kbinning.
num_bins; ibin++) {
646 powspec_out.
keff.push_back(k_save[ibin]);
647 powspec_out.
nmodes.push_back(nmodes_save[ibin]);
648 powspec_out.
pk_raw.push_back(norm_factor * pk_save[ibin]);
649 powspec_out.
pk_shot.push_back(norm_factor * sn_save[ibin]);
653 delete[] nmodes_save;
delete[] k_save;
delete[] pk_save;
delete[] sn_save;
657 "... computed power spectrum "
658 "from a periodic-box simulation-type catalogue "
659 "in the global plane-parallel approximation."
675 "Computing two-point correlation function "
676 "from a periodic-box simulation-type catalogue "
677 "in the global plane-parallel approximation."
685 int* npairs_save =
new int[rbinning.
num_bins];
686 double* r_save =
new double[rbinning.
num_bins];
687 std::complex<double>* xi_save =
new std::complex<double>[rbinning.
num_bins];
688 for (
int ibin = 0; ibin < rbinning.
num_bins; ibin++) {
689 npairs_save[ibin] = 0;
695 double norm = double(catalogue_data.
ntotal) * double(catalogue_data.
ntotal)
697 if (std::fabs(1 - norm * norm_factor) > eps_norm) {
700 "Two-point correlation function normalisation input differs from "
701 "expected value for an unweight field in a periodic box."
710#if defined(TRV_USE_OMP) && defined(TRV_USE_FFTWOMP)
719 std::complex<double> sn_amp = double(catalogue_data.
ntotal);
725 dn, dn, sn_amp, params.
ELL, 0, rbinning
728 for (
int ibin = 0; ibin < rbinning.
num_bins; ibin++) {
729 npairs_save[ibin] = stats_2pt.
npairs[ibin];
730 r_save[ibin] = stats_2pt.
r[ibin];
731 xi_save[ibin] += double(2*params.
ELL + 1) * stats_2pt.
xi[ibin];
740 for (
int ibin = 0; ibin < rbinning.
num_bins; ibin++) {
742 corrfunc_out.
reff.push_back(r_save[ibin]);
743 corrfunc_out.
npairs.push_back(npairs_save[ibin]);
744 corrfunc_out.
xi.push_back(norm_factor * xi_save[ibin]);
748 delete[] npairs_save;
delete[] r_save;
delete[] xi_save;
752 "... computed two-point correlation function "
753 "from a periodic-box simulation-type catalogue "
754 "in the global plane-parallel approximation."
764 double alpha,
double norm_factor
770 "Computing two-point correlation function window "
771 "from a random catalogue..."
780 int ell1 = params.
ELL;
783 int* npairs_save =
new int[rbinning.
num_bins];
784 double* r_save =
new double[rbinning.
num_bins];
785 std::complex<double>* xi_save =
new std::complex<double>[rbinning.
num_bins];
786 for (
int ibin = 0; ibin < rbinning.
num_bins; ibin++) {
787 npairs_save[ibin] = 0;
796#if defined(TRV_USE_OMP) && defined(TRV_USE_FFTWOMP)
800 MeshField dn_00(params,
true,
"`dn_00`");
806 for (
int M_ = - params.
ELL; M_ <= params.
ELL; M_++) {
807 MeshField dn_LM(params,
true,
"`dn_LM`");
809 catalogue_rand, los_rand, alpha, params.
ELL, M_
814 catalogue_rand, los_rand, alpha, params.
ELL, M_
819 for (
int m1 = - ell1; m1 <= ell1; m1++) {
824 dn_LM, dn_00, sn_amp, ell1, m1, rbinning
827 for (
int ibin = 0; ibin < rbinning.
num_bins; ibin++) {
828 xi_save[ibin] += coupling * stats_2pt.
xi[ibin];
831 if (M_ == 0 && m1 == 0) {
832 for (
int ibin = 0; ibin < rbinning.
num_bins; ibin++) {
833 npairs_save[ibin] = stats_2pt.
npairs[ibin];
834 r_save[ibin] = stats_2pt.
r[ibin];
841 "Two-point correlation function window term computed at order M = %d.",
852 for (
int ibin = 0; ibin < rbinning.
num_bins; ibin++) {
854 corrfunc_win_out.
reff.push_back(r_save[ibin]);
855 corrfunc_win_out.
npairs.push_back(npairs_save[ibin]);
856 corrfunc_win_out.
xi.push_back(norm_factor * xi_save[ibin]);
860 delete[] npairs_save;
delete[] r_save;
delete[] xi_save;
864 "... computed two-point correlation function window "
865 "from a random catalogue."
869 return corrfunc_win_out;
Isotropic coordinate binning.
std::vector< double > bin_centres
bin centres
int num_bins
number of bins
Field (pseudo-two-point) statistics.
std::vector< std::complex< double > > sn
shot-noise power in bins
std::vector< std::complex< double > > pk
pseudo power spectrum in bins
std::vector< int > nmodes
number of wavevector modes in bins
std::vector< double > k
average wavenumber in bins
std::vector< std::complex< double > > xi
pseudo two-point correlation function in bins
std::vector< int > npairs
number of separation pairs in bins
void compute_ylm_wgtd_2pt_stats_in_fourier(MeshField &field_a, MeshField &field_b, std::complex< double > shotnoise_amp, int ell, int m, trv::Binning &kbinning)
Compute binned two-point statistics in Fourier space.
void compute_ylm_wgtd_2pt_stats_in_config(MeshField &field_a, MeshField &field_b, std::complex< double > shotnoise_amp, int ell, int m, trv::Binning &rbinning)
Compute binned two-point statistics in configuration space.
Discretely sampled field on a mesh grid from particle catalogues.
void assign_weighted_field_to_mesh(ParticleCatalogue &particles, fftw_complex *weights)
Assign a weighted field to a mesh by interpolation scheme.
double calc_grid_based_powlaw_norm(ParticleCatalogue &particles, int order)
void fourier_transform()
Fourier transform the field.
void compute_ylm_wgtd_field(ParticleCatalogue &particles_data, ParticleCatalogue &particles_rand, LineOfSight *los_data, LineOfSight *los_rand, double alpha, int ell, int m)
Compute the weighted field (fluctuations) further weighted by the reduced spherical harmonics.
fftw_complex * field
complex field on mesh
void compute_unweighted_field_fluctuations_insitu(ParticleCatalogue &particles)
Compute the unweighted field fluctuations.
int validate()
Validate parameters.
long long nmesh
number of mesh grid cells
std::string assignment
mesh assignment scheme: {"ngp", "cic", "tsc" (default), "pcs"}
int ELL
spherical degree associated with the line of sight
double volume
box volume (in Mpc^3/h^3)
int ngrid[3]
grid number in each dimension
double boxsize[3]
box size (in Mpc/h) in each dimension
double pos_span[3]
span of particle coordinates
ParticleData * pdata
particle data
int ntotal
total number of particles
double wstotal
total sample weight of particles
Exception raised when the data to be operated on are invalid.
void error(const char *fmt_string,...)
Emit a warning-level message.
void warn(const char *fmt_string,...)
Emit a warning-level message.
void stat(const char *fmt_string,...)
Emit a status-level message.
void reset_level(LogLevel level)
Reset the logger threshold level.
const double eps_coupling
zero-tolerance for Wigner 3-j coupling coefficients
double wigner_3j(int j1, int j2, int j3, int m1, int m2, int m3)
Calculate Wigner 3-j symbol.
double gbytesMem
current memory usage in gibibytes
double size_in_gb(long long num)
Return size in gibibytes.
void update_maxmem()
Update the maximum memory usage estimate.
Logger logger
default logger (at NSET logging level)
double calc_powspec_shotnoise_from_particles(ParticleCatalogue &particles, double alpha=1.)
Calculate particle-based power spectrum shot noise level.
trv::TwoPCFWindowMeasurements compute_corrfunc_window(trv::ParticleCatalogue &catalogue_rand, trv::LineOfSight *los_rand, trv::ParameterSet ¶ms, trv::Binning rbinning, double alpha, double norm_factor)
Compute two-point correlation function window from a random catalogue and optionally save the results...
trv::TwoPCFMeasurements compute_corrfunc(ParticleCatalogue &catalogue_data, ParticleCatalogue &catalogue_rand, LineOfSight *los_data, LineOfSight *los_rand, trv::ParameterSet ¶ms, trv::Binning &rbinning, double norm_factor)
Compute two-point correlation function from paired survey-type catalogues.
double calc_coupling_coeff_2pt(int ell, int ELL, int m, int M)
Calculate the coupling coefficient for spherical-harmonic components of full two-point statistics.
std::complex< double > calc_ylm_wgtd_shotnoise_amp_for_powspec(ParticleCatalogue &particles_data, ParticleCatalogue &particles_rand, LineOfSight *los_data, LineOfSight *los_rand, double alpha, int ell, int m)
Calculate power spectrum shot noise weighted by reduced spherical harmonics.
trv::PowspecMeasurements compute_powspec_in_gpp_box(ParticleCatalogue &catalogue_data, trv::ParameterSet ¶ms, trv::Binning kbinning, double norm_factor)
Compute power spectrum in a periodic box in the global plane-parallel approximation.
double calc_powspec_normalisation_from_particles(ParticleCatalogue &particles, double alpha=1.)
Calculate particle-based power spectrum normalisation.
double calc_powspec_normalisation_from_mesh(trv::ParticleCatalogue &particles, trv::ParameterSet ¶ms, double alpha=1.)
Calculate mesh-based power spectrum normalisation.
trv::PowspecMeasurements compute_powspec(ParticleCatalogue &catalogue_data, ParticleCatalogue &catalogue_rand, LineOfSight *los_data, LineOfSight *los_rand, trv::ParameterSet ¶ms, trv::Binning &kbinning, double norm_factor)
Compute power spectrum from paired survey-type catalogues.
trv::TwoPCFMeasurements compute_corrfunc_in_gpp_box(ParticleCatalogue &catalogue_data, trv::ParameterSet ¶ms, trv::Binning &rbinning, double norm_factor)
Compute two-point correlation function in a periodic box in the global plane-parallel approximation.
double calc_powspec_normalisation_from_meshes(trv::ParticleCatalogue &particles_data, trv::ParticleCatalogue &particles_rand, trv::ParameterSet ¶ms, double alpha)
Calculate power spectrum normalisation from mixed meshes.
double pos[3]
3-d position vector
Power spectrum measurements.
int dim
dimension of data vector
std::vector< double > kbin
central wavenumber in bins
std::vector< int > nmodes
std::vector< std::complex< double > > pk_shot
power spectrum shot noise
std::vector< double > keff
effective wavenumber in bins
std::vector< std::complex< double > > pk_raw
power spectrum raw measurements (with normalisation and shot noise)
Two-point correlation function measurements.
std::vector< double > rbin
central separation in bins
std::vector< int > npairs
int dim
dimension of data vector
std::vector< std::complex< double > > xi
two-point correlation function measurements (with normalisation)
std::vector< double > reff
effective separation in bins
Two-point correlation function window measurements.
std::vector< std::complex< double > > xi
int dim
dimension of data vector
std::vector< double > reff
effective separation in bins
std::vector< double > rbin
central separation in bins
std::vector< int > npairs
Two-point statistic computations.