Triumvirate C++ API 0.5.0
Three-point clustering measurements in large-scale structure analyses.
Loading...
Searching...
No Matches
trv::MeshField Class Reference

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 &params, bool plan_ini=true, const std::string &name="mesh-field")
 Construct the mesh field.
 
 MeshField (trv::ParameterSet &params, 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
 

Detailed Description

Discretely sampled field on a mesh grid from particle catalogues.

Definition at line 68 of file field.hpp.

Constructor & Destructor Documentation

◆ MeshField() [1/2]

trv::MeshField::MeshField ( trv::ParameterSet & params,
bool plan_ini = true,
const std::string & name = "mesh-field" )
explicit

Construct the mesh field.

Parameters
paramsParameter set.
plan_iniFlag for FFTW plan initialisation (default is true).
nameField name (default is "mesh-field").

Definition at line 42 of file field.cpp.

+ Here is the call graph for this function:

◆ MeshField() [2/2]

trv::MeshField::MeshField ( trv::ParameterSet & params,
fftw_plan & transform,
fftw_plan & inv_transform,
const std::string & name = "mesh-field" )
explicit

Construct the mesh field with external FFTW plans.

Parameters
paramsParameter set.
transformExternal FFTW plan for Fourier transform.
inv_transformExternal FFTW plan for inverse Fourier transform.
nameField 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 273 of file field.cpp.

+ Here is the call graph for this function:

◆ ~MeshField()

trv::MeshField::~MeshField ( )

Destruct the mesh field.

Definition at line 329 of file field.cpp.

+ Here is the call graph for this function:

Member Function Documentation

◆ reset_density_field()

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 358 of file field.cpp.

+ Here is the caller graph for this function:

◆ operator[]()

const fftw_complex & trv::MeshField::operator[] ( long long gid)

Return mesh field grid cell value.

Parameters
gidGrid index.
Returns
Field value.

Definition at line 382 of file field.cpp.

◆ assign_weighted_field_to_mesh()

void trv::MeshField::assign_weighted_field_to_mesh ( ParticleCatalogue & particles,
fftw_complex * weights )

Assign a weighted field to a mesh by interpolation scheme.

Parameters
particlesParticle catalogue.
weightsWeight field.

Definition at line 427 of file field.cpp.

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ compute_unweighted_field()

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) \,. \]

Parameters
particlesParticle catalogue.

Definition at line 1066 of file field.cpp.

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ compute_unweighted_field_fluctuations_insitu()

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.

Parameters
particlesParticle catalogue.

Definition at line 1089 of file field.cpp.

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ compute_ylm_wgtd_field() [1/2]

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 \).

See also
Eq. (34) in Sugiyama et al. (2019) [1803.02132].
Parameters
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.
alphaAlpha contrast.
ellDegree of the spherical harmonic.
mOrder of the spherical harmonic.

Definition at line 1106 of file field.cpp.

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ compute_ylm_wgtd_field() [2/2]

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.

Parameters
particlesParticle catalogue.
losParticle lines of sight.
alphaAlpha contrast.
ellDegree of the spherical harmonic.
mOrder 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 1188 of file field.cpp.

+ Here is the call graph for this function:

◆ compute_ylm_wgtd_quad_field() [1/2]

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}} \,. \]

See also
Eq. (46) in Sugiyama et al. (2019) [1803.02132].
Parameters
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.
alphaAlpha contrast.
ellDegree of the spherical harmonic.
mOrder of the spherical harmonic.

Definition at line 1229 of file field.cpp.

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ compute_ylm_wgtd_quad_field() [2/2]

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.

Parameters
particlesParticle catalogue.
losParticle lines of sight.
alphaAlpha contrast.
ellDegree of the spherical harmonic.
mOrder 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 1316 of file field.cpp.

+ Here is the call graph for this function:

◆ fourier_transform()

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 1365 of file field.cpp.

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ inv_fourier_transform()

void trv::MeshField::inv_fourier_transform ( )

Inverse Fourier transform the (FFT-transformed) field.

Definition at line 1448 of file field.cpp.

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ apply_wide_angle_pow_law_kernel()

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 1479 of file field.cpp.

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ apply_assignment_compensation()

void trv::MeshField::apply_assignment_compensation ( )

Apply compensation needed after Fourier transform for assignment schemes.

Definition at line 1516 of file field.cpp.

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ inv_fourier_transform_ylm_wgtd_field_band_limited()

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}) \,. \]

See also
Eq. (42) in Sugiyama et al. (2019) [1803.02132].
Parameters
[in]field_fourierA Fourier-space field.
[in]ylmReduced spherical harmonic on a mesh.
[in]k_bandBand wavenumber.
[in]dk_bandBand wavenumber width.
[out]k_effEffective band wavenumber.
[out]nmodesNumber of wavevector modes in band.

Definition at line 1544 of file field.cpp.

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ inv_fourier_transform_sjl_ylm_wgtd_field()

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}) \,. \]

See also
Eq. (49) in Sugiyama et al. (2019) [1803.02132].
Parameters
field_fourierA Fourier-space field.
ylmReduced spherical harmonic on a mesh.
sjlSpherical Bessel function interpolator.
rSeparation in configuration space.

Definition at line 1627 of file field.cpp.

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ calc_grid_based_powlaw_norm()

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 \,. \]

Parameters
particles(Typically random-source) particle catalogue.
orderOrder \( N \) of the N-point statistics.
Returns
norm_factor Normalisation factor.

Definition at line 1702 of file field.cpp.

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

Friends And Related Symbol Documentation

◆ FieldStats

friend class FieldStats
friend

Definition at line 430 of file field.hpp.

Member Data Documentation

◆ params

trv::ParameterSet trv::MeshField::params

parameter set

Definition at line 70 of file field.hpp.

◆ name

std::string trv::MeshField::name

field name

Definition at line 71 of file field.hpp.

◆ field

fftw_complex* trv::MeshField::field

complex field on mesh

Definition at line 72 of file field.hpp.

◆ dr

double trv::MeshField::dr[3]

grid size in each dimension

Definition at line 73 of file field.hpp.

◆ dk

double trv::MeshField::dk[3]

fundamental wavenumber in each dimension

Definition at line 74 of file field.hpp.

◆ vol

double trv::MeshField::vol

mesh volume

Definition at line 75 of file field.hpp.

◆ vol_cell

double trv::MeshField::vol_cell

mesh grid cell volume

Definition at line 76 of file field.hpp.


The documentation for this class was generated from the following files: