|
Triumvirate C++ API 0.5.0.post1.dev301+g026f21751
Three-point clustering measurements in large-scale structure analyses.
|
Discretely sampled field on a mesh grid from particle catalogues. More...
#include <field.hpp>
Collaboration diagram for trv::MeshField:Public Member Functions | |
| MeshField (trv::ParameterSet ¶ms, bool plan_ini=true, const std::string &name="mesh-field") | |
| Construct the mesh field. | |
| MeshField (trv::ParameterSet ¶ms, fftw_plan &transform, fftw_plan &inv_transform, const std::string &name="mesh-field") | |
| Construct the mesh field with external FFTW plans. | |
| ~MeshField () | |
| Destruct the mesh field. | |
| void | reset_density_field () |
| (Re-)initialise the complex field (and its shadow) on mesh. | |
| const fftw_complex & | operator[] (long long gid) |
| Return mesh field grid cell value. | |
| void | assign_weighted_field_to_mesh (ParticleCatalogue &particles, fftw_complex *weights) |
| Assign a weighted field to a mesh by interpolation scheme. | |
| void | compute_unweighted_field (ParticleCatalogue &particles) |
| Compute the unweighted field. | |
| void | compute_unweighted_field_fluctuations_insitu (ParticleCatalogue &particles) |
| Compute the unweighted field fluctuations. | |
| 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. | |
| void | compute_ylm_wgtd_field (ParticleCatalogue &particles, LineOfSight *los, double alpha, int ell, int m) |
| Compute the weighted field further weighted by the reduced spherical harmonics. | |
| void | compute_ylm_wgtd_quad_field (ParticleCatalogue &particles_data, ParticleCatalogue &particles_rand, LineOfSight *los_data, LineOfSight *los_rand, double alpha, int ell, int m) |
| Compute the quadratic weighted field (fluctuations) further weighted by the reduced spherical harmonics. | |
| void | compute_ylm_wgtd_quad_field (ParticleCatalogue &particles, LineOfSight *los, double alpha, int ell, int m) |
| Compute the quadratic weighted field (fluctuations) further weighted by the reduced spherical harmonics. | |
| void | fourier_transform () |
| Fourier transform the field. | |
| void | inv_fourier_transform () |
| Inverse Fourier transform the (FFT-transformed) field. | |
| void | apply_wide_angle_pow_law_kernel () |
| Apply wide-angle correction kernel in configuration space. | |
| void | apply_assignment_compensation () |
| Apply compensation needed after Fourier transform for assignment schemes. | |
| void | inv_fourier_transform_ylm_wgtd_field_band_limited (MeshField &field_fourier, std::vector< std::complex< double > > &ylm, double k_band, double dk_band, double &k_eff, int &nmodes) |
| Inverse Fourier transform a field \( f \) weighted by the reduced spherical harmonics restricted to a wavenumber band. | |
| void | inv_fourier_transform_sjl_ylm_wgtd_field (MeshField &field_fourier, std::vector< std::complex< double > > &ylm, trvm::SphericalBesselCalculator &sjl, double r) |
| Inverse Fourier transform a field \( f \) weighted by the spherical Bessel function and reduced spherical harmonics. | |
| double | calc_grid_based_powlaw_norm (ParticleCatalogue &particles, int order) |
Public Attributes | |
| trv::ParameterSet | params |
| parameter set | |
| std::string | name |
| field name | |
| fftw_complex * | field |
| complex field on mesh | |
| double | dr [3] |
| grid size in each dimension | |
| double | dk [3] |
| fundamental wavenumber in each dimension | |
| double | vol |
| mesh volume | |
| double | vol_cell |
| mesh grid cell volume | |
Friends | |
| class | FieldStats |
Discretely sampled field on a mesh grid from particle catalogues.
|
explicit |
|
explicit |
Construct the mesh field with external FFTW plans.
| params | Parameter set. |
| transform | External FFTW plan for Fourier transform. |
| inv_transform | External FFTW plan for inverse Fourier transform. |
| name | Field name (default is "mesh-field"). |
This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
Definition at line 367 of file field.cpp.
Here is the call graph for this function:| trv::MeshField::~MeshField | ( | ) |
| void trv::MeshField::reset_density_field | ( | ) |
(Re-)initialise the complex field (and its shadow) on mesh.
This is an explicit method to reset values of trv::MeshField::field (and its interlaced counterpart) to zeros.
Definition at line 500 of file field.cpp.
Here is the caller graph for this function:| const fftw_complex & trv::MeshField::operator[] | ( | long long | gid | ) |
| void trv::MeshField::assign_weighted_field_to_mesh | ( | ParticleCatalogue & | particles, |
| fftw_complex * | weights ) |
| void trv::MeshField::compute_unweighted_field | ( | ParticleCatalogue & | particles | ) |
Compute the unweighted field.
This is is the number density field
\[ n(\vec{x}) = \sum_i \delta^{(\mathrm{D})}(\vec{x} - \vec{x}_i) \,. \]
| particles | Particle catalogue. |
Definition at line 1208 of file field.cpp.
Here is the call graph for this function:
Here is the caller graph for this function:| void trv::MeshField::compute_unweighted_field_fluctuations_insitu | ( | ParticleCatalogue & | particles | ) |
Compute the unweighted field fluctuations.
This is typically used for the number density fluctuations in a simulation snapshot in a periodic box, with a global field value calculated from the particle number and box volume subtracted to compute the fluctuations, i.e.
\[ \delta{n}(\vec{x}) = \sum_i \delta^{(\mathrm{D})}(\vec{x} - \vec{x}_i) - N / V \,, \]
where \( N \) is the total number of particles and \( V \) the box volume.
| particles | Particle catalogue. |
Definition at line 1229 of file field.cpp.
Here is the call graph for this function:
Here is the caller graph for this function:| void trv::MeshField::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.
For a field (or its fluctuations) \( f \), this is
\[ f_{LM}(\vec{x}) = {\sum_i}{\vphantom{\sum}}' y_{LM}^*(\hat{\vec{x}}) w(\vec{x}) \delta^{(\mathrm{D})}(\vec{x} - \vec{x}_i) \,, \]
where if a pair of catalogues are provided,
\[ {\sum_i}{\vphantom{\sum}}' = \sum_{i \in \mathrm{data}} - \alpha \sum_{i \in \mathrm{rand}} \]
with \( f \equiv \delta{n} \), and otherwise
\[ {\sum_i}{\vphantom{\sum}}' = \sum_{i \in \mathrm{data\ or\ rand}} \]
with \( f \equiv n \).
| particles_data | (Data-source) particle catalogue. |
| particles_rand | (Random-source) particle catalogue. |
| los_data | (Data-source) particle lines of sight. |
| los_rand | (Random-source) particle lines of sight. |
| alpha | Alpha contrast. |
| ell | Degree of the spherical harmonic. |
| m | Order of the spherical harmonic. |
Definition at line 1246 of file field.cpp.
Here is the call graph for this function:
Here is the caller graph for this function:| void trv::MeshField::compute_ylm_wgtd_field | ( | ParticleCatalogue & | particles, |
| LineOfSight * | los, | ||
| double | alpha, | ||
| int | ell, | ||
| int | m ) |
Compute the weighted field further weighted by the reduced spherical harmonics.
| particles | Particle catalogue. |
| los | Particle lines of sight. |
| alpha | Alpha contrast. |
| ell | Degree of the spherical harmonic. |
| m | Order of the spherical harmonic. |
This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
Definition at line 1325 of file field.cpp.
Here is the call graph for this function:| void trv::MeshField::compute_ylm_wgtd_quad_field | ( | ParticleCatalogue & | particles_data, |
| ParticleCatalogue & | particles_rand, | ||
| LineOfSight * | los_data, | ||
| LineOfSight * | los_rand, | ||
| double | alpha, | ||
| int | ell, | ||
| int | m ) |
Compute the quadratic weighted field (fluctuations) further weighted by the reduced spherical harmonics.
For a quadratic field (of fluctuations) \( f \), this is
\[ f_{LM}(\vec{x}) = {\sum_i}{\vphantom{\sum}}' y_{LM}^*(\hat{\vec{x}}) w(\vec{x})^2 \delta^{(\mathrm{D})}(\vec{x} - \vec{x}_i) \,, \]
where if a pair of catalogues are provided,
\[ {\sum_i}{\vphantom{\sum}}' = \sum_{i \in \mathrm{data}} + \alpha^2 \sum_{i \in \mathrm{rand}} \,, \]
and otherwise
\[ {\sum_i}{\vphantom{\sum}}' = \sum_{i \in \mathrm{data\ or\ rand}} \,. \]
| particles_data | (Data-source) particle catalogue. |
| particles_rand | (Random-source) particle catalogue. |
| los_data | (Data-source) particle lines of sight. |
| los_rand | (Random-source) particle lines of sight. |
| alpha | Alpha contrast. |
| ell | Degree of the spherical harmonic. |
| m | Order of the spherical harmonic. |
Definition at line 1364 of file field.cpp.
Here is the call graph for this function:
Here is the caller graph for this function:| void trv::MeshField::compute_ylm_wgtd_quad_field | ( | ParticleCatalogue & | particles, |
| LineOfSight * | los, | ||
| double | alpha, | ||
| int | ell, | ||
| int | m ) |
Compute the quadratic weighted field (fluctuations) further weighted by the reduced spherical harmonics.
| particles | Particle catalogue. |
| los | Particle lines of sight. |
| alpha | Alpha contrast. |
| ell | Degree of the spherical harmonic. |
| m | Order of the spherical harmonic. |
This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
Definition at line 1449 of file field.cpp.
Here is the call graph for this function:| void trv::MeshField::fourier_transform | ( | ) |
Fourier transform the field.
If trv::MeshField.params.interlace is set to "true", interlacing is performed where a phase factor is multiplied into the 'shadow' complex field before the average of the complex field and its shadow is taken.
Definition at line 1496 of file field.cpp.
Here is the call graph for this function:
Here is the caller graph for this function:| void trv::MeshField::inv_fourier_transform | ( | ) |
| void trv::MeshField::apply_wide_angle_pow_law_kernel | ( | ) |
Apply wide-angle correction kernel in configuration space.
This multiplies the field by a power law \( r^{- i - j} \) at order \( (i, j) \) where \( r \) is the grid cell radial distance.
Definition at line 1727 of file field.cpp.
Here is the call graph for this function:
Here is the caller graph for this function:| void trv::MeshField::apply_assignment_compensation | ( | ) |
| void trv::MeshField::inv_fourier_transform_ylm_wgtd_field_band_limited | ( | MeshField & | field_fourier, |
| std::vector< std::complex< double > > & | ylm, | ||
| double | k_band, | ||
| double | dk_band, | ||
| double & | k_eff, | ||
| int & | nmodes ) |
Inverse Fourier transform a field \( f \) weighted by the reduced spherical harmonics restricted to a wavenumber band.
This method computes the quantity
\[ F_{LM}(\vec{x}; k) = \frac{(2\pi)^3}{4\pi k^2} \int \frac{\mathrm{d}^3\,\vec{k}'}{(2\pi)^3} \mathrm{e}^{\mathrm{i} \vec{k}' \cdot \vec{x}} \delta^{(\mathrm{D})}(k' - k) y_{LM}(\hat{\vec{k}}) f(\vec{k}) \,. \]
| [in] | field_fourier | A Fourier-space field. |
| [in] | ylm | Reduced spherical harmonic on a mesh. |
| [in] | k_band | Band wavenumber. |
| [in] | dk_band | Band wavenumber width. |
| [out] | k_eff | Effective band wavenumber. |
| [out] | nmodes | Number of wavevector modes in band. |
Definition at line 1792 of file field.cpp.
Here is the call graph for this function:
Here is the caller graph for this function:| void trv::MeshField::inv_fourier_transform_sjl_ylm_wgtd_field | ( | MeshField & | field_fourier, |
| std::vector< std::complex< double > > & | ylm, | ||
| trvm::SphericalBesselCalculator & | sjl, | ||
| double | r ) |
Inverse Fourier transform a field \( f \) weighted by the spherical Bessel function and reduced spherical harmonics.
This method computes the quantity
\[ F_{LM}(\vec{x}; r) = \mathrm{i}^L \int \frac{\mathrm{d}^3\,\vec{k}}{(2\pi)^3} \mathrm{e}^{\mathrm{i} \vec{k} \cdot \vec{x}} j_L(kr) y_{LM}(\hat{\vec{k}}) f(\vec{k}) \,. \]
| field_fourier | A Fourier-space field. |
| ylm | Reduced spherical harmonic on a mesh. |
| sjl | Spherical Bessel function interpolator. |
| r | Separation in configuration space. |
Definition at line 1908 of file field.cpp.
Here is the call graph for this function:
Here is the caller graph for this function:| double trv::MeshField::calc_grid_based_powlaw_norm | ( | ParticleCatalogue & | particles, |
| int | order ) |
Calculate the normalisation factor \( 1/I_N \) for N-point (for two-point statistics), where
\[ I_N = \int\mathrm{d}^3\,\vec{x} w(\vec{x})^N \bar{n}(\vec{x})^N \,. \]
| particles | (Typically random-source) particle catalogue. |
| order | Order \( N \) of the N-point statistics. |
Definition at line 2017 of file field.cpp.
Here is the call graph for this function:
Here is the caller graph for this function:
|
friend |
| trv::ParameterSet trv::MeshField::params |
| double trv::MeshField::dk[3] |