CosmoBolognaLib
Free Software C++/Python libraries for cosmological calculations
cbl::measure::angularpk::PowerSpectrum_angular Class Reference

The class PowerSpectrum_angular. More...

#include <PowerSpectrum_Angular.h>

Inheritance diagram for cbl::measure::angularpk::PowerSpectrum_angular:
Collaboration diagram for cbl::measure::angularpk::PowerSpectrum_angular:

Public Member Functions

void set_l_output (const BinType binType)
 generate vector which contains the multipoles at which the angular power spectrum is computed More...
 
void measure (const AngularEstimator estimator=AngularEstimator::_Fast_, const std::string dir_correlation_input="", const std::string file_correlation_input="", const int n_lines_header=1, CoordinateUnits inputUnits=CoordinateUnits::_arcminutes_, BinType input_binType=BinType::_logarithmic_, cbl::measure::ErrorType errorType=ErrorType::_Poisson_, const std::string dir_correlation_output=par::defaultString, const std::string file_correlation_output="xi_angular.dat", const int nMocks=0)
 measure the angular power spectrum More...
 
void measureFast (const std::string dir_correlation_input="", const std::string file_correlation_input="", const int n_lines_header=1, CoordinateUnits inputUnits=CoordinateUnits::_arcminutes_, BinType input_binType=BinType::_logarithmic_, cbl::measure::ErrorType errorType=ErrorType::_Poisson_, const std::string dir_correlation_output=par::defaultString, const std::string file_correlation_output="xi_angular.dat", const int nMocks=0)
 measure the angular power spectrum with fast estimator More...
 
void measureSphericalArmonic ()
 measure the angular power spectrum with spherical armonic estimator
 
std::vector< double > theta ()
 get the private member m_theta More...
 
std::vector< double > Wl ()
 get the private member Wl More...
 
std::vector< double > AngularPowerSpectrum ()
 get the private member AngularPowerSpectrum More...
 
std::vector< double > l ()
 get the private member m_l More...
 
void compute_mixing_matrix (std::string dir_window_input="", std::string file_window_input="")
 compute the mixing matrix More...
 
double angular_mask (double theta, double RA)
 include the binary angular mask More...
 
Constructors/destructors
 PowerSpectrum_angular ()=default
 default constructor Power_spectrum_angular
 
 PowerSpectrum_angular (const catalogue::Catalogue data, const catalogue::Catalogue random, const double l_min, const double l_max, const int Nl, const BinType binType, const double thetaMin, const double thetaMax, const int nbins, const double shift, const CoordinateUnits angularUnits, const BinType correlation_binType)
 constructor More...
 
 PowerSpectrum_angular (const catalogue::Catalogue data, const double l_min, const double l_max, const int bandwidth=1, const std::string mask_file="", const std::string mask_type="", const double pixel_area=0, const int n_lines_header=1)
 constructor More...
 
 ~PowerSpectrum_angular ()=default
 default destructor
 
Input/Output methods
void write (const std::string dir, const std::string file)
 write the angular power spectrum on file More...
 
void write_average (const std::string dir, const std::string file)
 write the angular power spectrum on file More...
 
void write_mixing_matrix (const std::string dir, const std::string file, bool store_window=false)
 write the angular power spectrum on file More...
 
void read (const std::string dir, const std::string file, const std::vector< int > column_x={1}, const std::vector< int > column_y={2}, const std::vector< int > column_error={3}, const int n_lines_header=1)
 read data from file More...
 
- Public Member Functions inherited from cbl::measure::Measure
 Measure ()=default
 default constructor
 
virtual ~Measure ()=default
 default destructor
 
virtual std::shared_ptr< data::Datadataset () const
 get the protected member dataset More...
 

Private Member Functions

void m_convert_angular_units (cbl::CoordinateUnits inputUnits)
 converts the input angle in radians More...
 
double m_dtheta_theta (const BinType binType, int i)
 compute the dtheta*theta for the integral More...
 
std::vector< double > m_error (const BinType binType, std::vector< double > w_err)
 measure the angular power spectrum errors More...
 
void set_catalogue (cbl::catalogue::Catalogue catalogue)
 set the catalogue for Spherical armonic estimator More...
 
void m_averagePowerSpectrum ()
 compute the average of the power spectrum
 
void m_read_angular_mask (const std::string mask_file, std::string mask_type, double pixel_area, int n_lines_header=1)
 read the angular mask from file More...
 

Private Attributes

std::vector< double > m_theta
 vector containing the angular separation
 
std::vector< double > m_w
 vector containing the angular correlation function
 
std::vector< double > m_w_err
 vector containing the angular correlation function error
 
std::vector< double > m_x
 vector containing the x coordinate read from file
 
std::vector< double > m_y
 vector containing the y coordinate read from file
 
std::vector< double > m_y_err
 vector containing the y errors read from file
 
std::vector< double > m_l
 vector containing the multipoles at which power spectrum is computed
 
std::vector< double > m_AngularPowerSpectrum
 vector containing the angular power spectrum
 
std::vector< double > m_l_av
 vector containing the average of the multipoles
 
std::vector< double > m_AngularPowerSpectrum_av
 vector containing the angular power spectrum average
 
double m_l_min
 the minimum multipole
 
double m_l_max
 the maximum multipole
 
double m_Nl
 the number of multipoles at which the angular power spectrum is computed
 
double m_nbins
 the number of multipoles at which the angular power spectrum is computed
 
std::vector< double > m_Wl
 vector containing the angular power spectrum of the window function
 
std::vector< std::vector< double > > m_RR
 the mixing matrix
 
std::string m_dir_correlation_input
 the directory to read w
 
std::string m_file_correlation_input
 the file to read w
 
int m_n_lines_header
 the header lines to skip
 
cbl::CoordinateUnits m_inputUnits
 the coordinate units of the input file
 
std::string m_dir_correlation_output
 the directory to write w (if computed with cbl)
 
std::string m_file_correlation_output
 the file to write w (if computed with cbl)
 
cbl::CoordinateUnits m_angularUnits
 the coordinate units of the correlation function (if computed with cbl)
 
cbl::BinType m_binType
 the bin type of the correlation function (if computed with cbl)
 
std::shared_ptr< twopt::TwoPointCorrelation1D_angularm_TwoPointCorrelation1D_angular
 variable that store the constructor for the measure of angular correlation function
 
std::shared_ptr< catalogue::Cataloguem_catalogue
 variable that store the constructor for the catalogue
 
std::vector< double > m_theta_mask
 vector containing the colatitude vector of the mask
 
std::vector< double > m_theta_mask_min
 vector containing the min colatitude vector of the mask
 
std::vector< double > m_theta_mask_max
 vector containing the max colatitude vector of the mask
 
std::vector< double > m_RA_mask
 vector containing the RA vector of the mask
 
std::vector< double > m_RA_mask_min
 vector containing the min RA vector of the mask
 
std::vector< double > m_RA_mask_max
 vector containing the max RA vector of the mask
 
double m_pixel_area
 pixel area
 
double m_survey_area
 the survey area
 

Additional Inherited Members

- Protected Attributes inherited from cbl::measure::Measure
std::shared_ptr< data::Datam_dataset
 the dataset of the measure
 

Detailed Description

The class PowerSpectrum_angular.

This class is used to handle objects of type PowerSpectrum_angular . It is used to measure the angular power spectrum, \(C_l\)

Examples
power_spectrum_angular.cpp.

Definition at line 89 of file PowerSpectrum_Angular.h.

Constructor & Destructor Documentation

◆ PowerSpectrum_angular() [1/2]

cbl::measure::angularpk::PowerSpectrum_angular::PowerSpectrum_angular ( const catalogue::Catalogue  data,
const catalogue::Catalogue  random,
const double  l_min,
const double  l_max,
const int  Nl,
const BinType  binType,
const double  thetaMin,
const double  thetaMax,
const int  nbins,
const double  shift,
const CoordinateUnits  angularUnits,
const BinType  correlation_binType 
)

constructor

Parameters
dataobject of class Catalogue containing the input catalogue
randomof class Catalogue containing the random data catalogue
l_minMinimum angular power spectrum multipole
l_maxMaximum angular power spectrum multipole
Nlnumber of multipoles
binTypebinning type
thetaMinminimum angular separation used to count the pairs
thetaMaxmaximum angular separation used to count the pairs
nbinsnumber of bins
shiftshift parameter, i.e. the radial shift is binSize*shift
angularUnitsangular units (of the output angular correlation file)
correlation_binTypethe bin type of the correlation function

Definition at line 176 of file AngularPowerSpectrum.cpp.

◆ PowerSpectrum_angular() [2/2]

cbl::measure::angularpk::PowerSpectrum_angular::PowerSpectrum_angular ( const catalogue::Catalogue  data,
const double  l_min,
const double  l_max,
const int  bandwidth = 1,
const std::string  mask_file = "",
const std::string  mask_type = "",
const double  pixel_area = 0,
const int  n_lines_header = 1 
)

constructor

Parameters
dataobject of class Catalogue containing the input catalogue
l_minMinimum angular power spectrum multipole
l_maxMaximum angular power spectrum multipole
bandwidththe bandwidth
mask_filethe file containing the pixel center positions of the mask in colatitude-RA for every pixel (obtained with Healpix pixelfunc.pix2ang), or theta, theta_min, theta_max, ra, ra_min, ra_max and pixel_area, for every pixel, i.e. pixel centers, borders and area. Leave it empty if there is not a mask file, in this case the survey area will be computed as \( (cos(\theta_{min})-cos(\theta_{max}))*(RA_{max}-RA_{min})\)
mask_typethe mask_type. "Healpix" if the mask file is obtained with the healpix function pixelfunc.pix2ang. Anything else if themask is obtained in other way
pixel_areathe average pixel area
n_lines_headerthe number of lines to skip when reading the mask file

Definition at line 195 of file AngularPowerSpectrum.cpp.

Member Function Documentation

◆ angular_mask()

double cbl::measure::angularpk::PowerSpectrum_angular::angular_mask ( double  theta,
double  RA 
)

include the binary angular mask

Parameters
thetathe colatitude to evaluate
RAthe right ascension to evaluate
Returns
0 if theta, RA is within a masked pixel, 1 otherwise

Definition at line 145 of file AngularPowerSpectrum.cpp.

◆ AngularPowerSpectrum()

std::vector<double> cbl::measure::angularpk::PowerSpectrum_angular::AngularPowerSpectrum ( )
inline

get the private member AngularPowerSpectrum

Returns
the angular power spectrum vector

Definition at line 439 of file PowerSpectrum_Angular.h.

◆ compute_mixing_matrix()

void cbl::measure::angularpk::PowerSpectrum_angular::compute_mixing_matrix ( std::string  dir_window_input = "",
std::string  file_window_input = "" 
)

compute the mixing matrix

Parameters
dir_window_inputinput directory harmonic coefficients of the mask
file_window_inputinput file harmonic coefficients of the mask

Definition at line 622 of file AngularPowerSpectrum.cpp.

◆ l()

std::vector<double> cbl::measure::angularpk::PowerSpectrum_angular::l ( )
inline

get the private member m_l

Returns
the multipoles vector

Definition at line 445 of file PowerSpectrum_Angular.h.

◆ m_convert_angular_units()

void cbl::measure::angularpk::PowerSpectrum_angular::m_convert_angular_units ( cbl::CoordinateUnits  inputUnits)
private

converts the input angle in radians

Parameters
inputUnitsinput units

Definition at line 425 of file AngularPowerSpectrum.cpp.

◆ m_dtheta_theta()

double cbl::measure::angularpk::PowerSpectrum_angular::m_dtheta_theta ( const BinType  binType,
int  i 
)
private

compute the dtheta*theta for the integral

Parameters
binTypethe binning type, linear or logarithmic
iindex to select the angle theta[i]
Returns
the dtheta*theta product

Definition at line 227 of file AngularPowerSpectrum.cpp.

◆ m_error()

std::vector< double > cbl::measure::angularpk::PowerSpectrum_angular::m_error ( const BinType  binType,
std::vector< double >  w_err 
)
private

measure the angular power spectrum errors

Parameters
binTypethe binning type, linear or logarithmic
w_errangular correlation function vector error
Returns
the error vector

Definition at line 400 of file AngularPowerSpectrum.cpp.

◆ m_read_angular_mask()

void cbl::measure::angularpk::PowerSpectrum_angular::m_read_angular_mask ( const std::string  mask_file,
std::string  mask_type,
double  pixel_area,
int  n_lines_header = 1 
)
private

read the angular mask from file

Parameters
mask_filethe file with pixel mask positions
mask_type"Healpix" in case of binary mask produced with pixelfunc.pix2ang in healpix. Else other case
pixel_areathe pixel area of healpix mask
n_lines_headerthe header lines to skip

Definition at line 51 of file AngularPowerSpectrum.cpp.

◆ measure()

void cbl::measure::angularpk::PowerSpectrum_angular::measure ( const AngularEstimator  estimator = AngularEstimator::_Fast_,
const std::string  dir_correlation_input = "",
const std::string  file_correlation_input = "",
const int  n_lines_header = 1,
CoordinateUnits  inputUnits = CoordinateUnits::_arcminutes_,
BinType  input_binType = BinType::_logarithmic_,
cbl::measure::ErrorType  errorType = ErrorType::_Poisson_,
const std::string  dir_correlation_output = par::defaultString,
const std::string  file_correlation_output = "xi_angular.dat",
const int  nMocks = 0 
)

measure the angular power spectrum

Parameters
estimatorthe power spectrum estimator
dir_correlation_inputAngular correlation function input directory
file_correlation_inputAngular correlation function input file
n_lines_headerthe header lines to skip
inputUnitsangular units in input file
input_binTypethe binning type, linear or logarithmic
errorTypethe error type, Poisson, Jackknife or Bootstrap
dir_correlation_outputAngular correlation function ouput directory (if measured by CBL)
file_correlation_outputAngular correlation function ouput file (if measured by CBL)
nMocksnumber of mocks for bootstrap
Examples
power_spectrum_angular.cpp.

Definition at line 451 of file AngularPowerSpectrum.cpp.

◆ measureFast()

void cbl::measure::angularpk::PowerSpectrum_angular::measureFast ( const std::string  dir_correlation_input = "",
const std::string  file_correlation_input = "",
const int  n_lines_header = 1,
CoordinateUnits  inputUnits = CoordinateUnits::_arcminutes_,
BinType  input_binType = BinType::_logarithmic_,
cbl::measure::ErrorType  errorType = ErrorType::_Poisson_,
const std::string  dir_correlation_output = par::defaultString,
const std::string  file_correlation_output = "xi_angular.dat",
const int  nMocks = 0 
)

measure the angular power spectrum with fast estimator

Parameters
dir_correlation_inputAngular correlation function input directory
file_correlation_inputAngular correlation function input file
n_lines_headerthe header lines to skip
inputUnitsangular units in input file
input_binTypethe binning type, linear or logarithmic
errorTypethe error type, Poisson, Jackknife or Bootstrap
dir_correlation_outputAngular correlation function ouput directory (if measured by CBL)
file_correlation_outputAngular correlation function ouput file (if measured by CBL)
nMocksnumber of mocks for bootstrap

Definition at line 471 of file AngularPowerSpectrum.cpp.

◆ read()

void cbl::measure::angularpk::PowerSpectrum_angular::read ( const std::string  dir,
const std::string  file,
const std::vector< int >  column_x = {1},
const std::vector< int >  column_y = {2},
const std::vector< int >  column_error = {3},
const int  n_lines_header = 1 
)

read data from file

Parameters
diroutput directory
fileoutput file
column_xvector of x column
column_yvector of y column
column_errorvector of error column
n_lines_headerthe header lines to skip

Definition at line 352 of file AngularPowerSpectrum.cpp.

◆ set_catalogue()

void cbl::measure::angularpk::PowerSpectrum_angular::set_catalogue ( cbl::catalogue::Catalogue  catalogue)
private

set the catalogue for Spherical armonic estimator

Parameters
cataloguethe catalogue to set

Definition at line 216 of file AngularPowerSpectrum.cpp.

◆ set_l_output()

void cbl::measure::angularpk::PowerSpectrum_angular::set_l_output ( const BinType  binType)

generate vector which contains the multipoles at which the angular power spectrum is computed

Parameters
binTypethe binning type, linear or logarithmic

Definition at line 243 of file AngularPowerSpectrum.cpp.

◆ theta()

std::vector<double> cbl::measure::angularpk::PowerSpectrum_angular::theta ( )
inline

get the private member m_theta

Returns
the angular separation vector

Definition at line 427 of file PowerSpectrum_Angular.h.

◆ Wl()

std::vector<double> cbl::measure::angularpk::PowerSpectrum_angular::Wl ( )
inline

get the private member Wl

Returns
the angular power spectrum of the window function vector

Definition at line 433 of file PowerSpectrum_Angular.h.

◆ write()

void cbl::measure::angularpk::PowerSpectrum_angular::write ( const std::string  dir,
const std::string  file 
)

write the angular power spectrum on file

Parameters
diroutput directory
fileoutput file

Definition at line 269 of file AngularPowerSpectrum.cpp.

◆ write_average()

void cbl::measure::angularpk::PowerSpectrum_angular::write_average ( const std::string  dir,
const std::string  file 
)

write the angular power spectrum on file

Parameters
diroutput directory
fileoutput file

Definition at line 293 of file AngularPowerSpectrum.cpp.

◆ write_mixing_matrix()

void cbl::measure::angularpk::PowerSpectrum_angular::write_mixing_matrix ( const std::string  dir,
const std::string  file,
bool  store_window = false 
)

write the angular power spectrum on file

Parameters
diroutput directory
fileoutput file
store_windowtrue → store the harmonic coefficients of the mask; false → do not store the harmonic coefficients of the mask

Definition at line 319 of file AngularPowerSpectrum.cpp.


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