API

The main class for handling the power measurements is the PowerData class. To initialize this class, the parameters are:

covariance The string specifying the name of the file holding the covariance matrix.
covariance_Nmocks The number of mocks that was used to measure the covariance matrix.
covariance_rescaling Rescale the covariance matrix read from file by this amount.
data_file The string specifying the name of the file holding the data measurements.
ells A list of integers specifying multipole numbers for each statistic in the final analysis.
fitting_range The \(k\) fitting range for each statistics.
grid_file A string specifying the name of the file holding a pyRSD.rsd.transfers.PkmuGrid to read.
max_ellprime When convolving a multipole of order ell, include contributions up to and including this number.
mode The type of data, either pkmu or poles
mu_bounds A list of tuples specifying the edges of the \(\mu\) bins.
statistics A list of the string names for each statistic that will be read from file
usedata A list of the statistic numbers that will be included in the final analysis.
window_file A string specifying the name of the file holding the correlation function multipoles of the window function.
class pyRSD.rsdfit.data.PowerData(param_file)

Class to hold several PowerMeasurement objects and combine the associated covariance matrices

Note

See the PowerData.help() function for a list of the parameters needed to initialize this class

Parameters:

param_file : str

the name of the parameter file holding the necessary parameters needed to initialize the object

covariance

The string specifying the name of the file holding the covariance matrix.

covariance_Nmocks

The number of mocks that was used to measure the covariance matrix.

If this is non-zero, then the inverse covariance matrix will be rescaled to account for noise due to the finite number of mocks

covariance_rescaling

Rescale the covariance matrix read from file by this amount.

data_file

The string specifying the name of the file holding the data measurements.

ells

A list of integers specifying multipole numbers for each statistic in the final analysis.

This must be supplied when the mode is poles

fitting_range

The \(k\) fitting range for each statistics.

This can either be a tuple of (kmin, kmax), which will be used for each statistic or a list of tuples of (kmin, kmax)

grid_file

A string specifying the name of the file holding a pyRSD.rsd.transfers.PkmuGrid to read.

static help()

Print out the help information for the necessary initialization parameters

max_ellprime

When convolving a multipole of order ell, include contributions up to and including this number.

mode

The type of data, either pkmu or poles

mu_bounds

A list of tuples specifying the edges of the \(\mu\) bins.

This should have (mu_min, mu_max), corresponding to the edges of the bins for each statistic in the final analysis

This must be supplied when the mode is pkmu

statistics

A list of the string names for each statistic that will be read from file

These strings should be of the form:

>> [‘pole_0’, ‘pole_2’, …] >> [‘pkmu_0.1’, ‘pkmu_0.3’, …]

to_file(filename, mode='w')

Save the parameters of this data class to a file

Parameters:

filename : str

the name of the file to write out the parameters to

mode : str

the mode to use when writing the parameters to file; i.e., ‘w’, ‘a’

usedata

A list of the statistic numbers that will be included in the final analysis.

This allows the user to exclude certain statistics read from file. By default (None), all statistics are included

window_file

A string specifying the name of the file holding the correlation function multipoles of the window function.

The file should contain columns of data, with the first column specifying the separation array \(s\), and the other columns giving the even-numbered correlation function multipoles of the window

Power Statistics

class pyRSD.rsdfit.data.PowerMeasurements

A list of PowerMeasurement objects

classmethod from_array(names, data)

Create a PowerMeasurement from the input structured array of data

Parameters:

names : list of str

the list of names for each measurement to load; each name must begin with pkmu_ or pole_ depending on the type of data. Examples of names include pkmu_0.1 for a \(\mu=0.1\) \(P(k,\mu)\) bin, or pole_0 for the monopole, pole_2 for the quadrupole, etc

data : array_like

a structured array holding the data fields; this have k, mu, power fields for \(P(k,\mu)\) data or k, power fields for multipole data

classmethod from_plaintext(names, filename)

Load a set of power measurements from a plaintex file

Parameters:

names : list of str

the list of names for each measurement to load; each name must begin with pkmu_ or pole_ depending on the type of data. Examples of names include pkmu_0.1 for a \(\mu=0.1\) \(P(k,\mu)\) bin, or pole_0 for the monopole, pole_2 for the quadrupole, etc

filename : str

the name of the file to load to load a structured array of data from

to_plaintext(filename)

Write out the data from the power measurements to a plaintext file

Covariance Matrix

class pyRSD.rsdfit.data.PoleCovarianceMatrix(data, k_center, ells, **kwargs)

Class to hold a covariance matrix for multipole measurements

Parameters:

data : array_like (N, ) or (N, N)

The data representing the covariance matrix. If 1D, the input is interpreted as the diagonal elements, with all other elements zero

k_center : array_like, (N, )

The wavenumbers where the center of the measurement k-bins are defined, for each element of the data matrix along one axis

ells : array_like, (N, )

The multipole numbers for each element of the data matrix along one axis

Examples

>>> from pyRSD.rsdfit.data import PoleCovarianceMatrix
>>> import numpy as np
>>> # generate 20 k bins from 0.01 to 0.4
>>> Nk = 20
>>> k_edges = np.linspace(0.01, 0.4, Nk+1)
>>> k_cen = 0.5 * (k_edges[1:] + k_edges[:-1])
>>> # generate 100 fake monopoles
>>> P0 = np.random.random(size=(100, Nk)) # Nk = 20, see above
>>> # 100 fake quadrupoles
>>> P2 = np.random.random(size=(100, Nk))
>>> # 100 fake hexadecapoles
>>> P4 = np.random.random(size=(100, Nk))
>>> # make the full data vector
>>> D = np.concatenate([P0, P2, P4], axis=-1) # shape is (100, 60)
>>> # compute the covariance matrix
>>> cov = np.cov(D, rowvar=False) # shape is (20,20)
>>> # initialize the PoleCovarianceMatrix
>>> ells = [0,2,4]
>>> C = PoleCovarianceMatrix(cov, k_cen, ells)
classmethod cutsky_gaussian_covariance(model, k, ells, nbar, fsky, zmin, zmax, FKP_P0=10000.0, Nmu=100, Nz=50)

Return the Gaussian prediction for the covariance between multipoles for a “cutsky” survey, i.e., a survey with a varying \(n(z)\) distribution

Parameters:

model : GalaxySpectrum, QuasarSpectrum

the model instance used to evaluate the theoretical \(P(k,\mu)\)

k : array_like

the array of wavenumbers (units of \(h/\mathrm{Mpc}\)) where the covariance matrix will be evaluated

ells : list of int

the list of multipole numbers to compute the covariance of

nbar : callable

a callable function returning the number density as a function of redshift

fsky : float

the sky fraction, used to normalized the effective volume calculation

zmin : float

the minimum redshift value to include when performing the effective volume calculation

zmax : float

the maximum redshift value to include when performing the effective volume calculation

FKP_P0 : float, optional

the FKP P0 value to use, where the FKP weights are \(1/(1+n(z)P_0)\); default is 1e4

Nmu : int, optional

the number of mu bins to use when performing the multipole integration over \(\mu\)

Nz : int, optional

the number of redshift bins to use when performing the integral over redshift

Returns:

PoleCovarianceMatrix :

the covariance matrix object holding the Gaussian prediction for the covariance between the specified multipoles

Examples

>>> import numpy
>>> from scipy.interpolate import InterpolatedUnivariateSpline as spline
>>> from pyRSD.rsd import GalaxySpectrum
>>> from pyRSD.rsdfit.data import PoleCovarianceMatrix
>>> filename = 'nbar_DR12v5_CMASSLOWZ_North_om0p31_Pfkp10000.dat'
>>> nbar = numpy.loadtxt(filename, skiprows=3)
>>> nbar = spline(nbar[:,0], nbar[:,3])
>>> fsky = 0.1436
>>> model = GalaxySpectrum(params='boss_dr12_fidcosmo.ini')
>>> k = numpy.arange(0., 0.4, 0.005) + 0.005/2
>>> ells = [0,2,4]
>>> zmin = 0.2
>>> zmax = 0.5
>>> C = PoleCovarianceMatrix.cutsky_gaussian_covariance(model, k, ells, nbar, fsky, zmin, zmax)
classmethod from_plaintext(filename)

Load the covariance matrix from a plain text file

classmethod periodic_gaussian_covariance(model, k, ells, nbar, volume, Nmu=100)

Return the Gaussian prediction for the covariance between multipoles for a periodic box simulation, where the number density is constant.

See eq. 16 of Grieb et al. 2015 arxiv:1509.04293

Parameters:

model : GalaxySpectrum, QuasarSpectrum

the model instance used to evaluate the theoretical \(P(k,\mu)\)

k : array_like

the array of wavenumbers (units of \(h/\mathrm{Mpc}\)) where the covariance matrix will be evaluated

ells : list of int

the list of multipole numbers to compute the covariance of

nbar : float

the constant number density in the box in units of \((\mathrm{Mpc}/h)^{-3}\) – the shot noise contribution to the covariance is the inverse of this value

volume : float

the volume of the box, in units of \((\mathrm{Mpc}/h)^3\)

Nmu : int, optional

the number of mu bins to use when performing the multipole integration over \(\mu\)

Returns:

PoleCovarianceMatrix :

the covariance matrix object holding the Gaussian prediction for the covariance between the specified multipoles

Examples

>>> import numpy
>>> from pyRSD.rsd import GalaxySpectrum
>>> from pyRSD.rsdfit.data import PoleCovarianceMatrix
>>> volume = 1380.0**3
>>> nbar = 3e-4
>>> model = GalaxySpectrum(params='boss_dr12_fidcosmo.ini')
>>> k = numpy.arange(0., 0.4, 0.005) + 0.005/2
>>> ells = [0,2,4]
>>> C = PoleCovarianceMatrix.periodic_gaussian_covariance(model, k, ells, nbar, volume)
to_plaintext(filename)

Save the covariance matrix to a plain text file

class pyRSD.rsdfit.data.PkmuCovarianceMatrix(data, k_center, mu_center, **kwargs)

Class to hold a covariance matrix for P(k, mu)

Parameters:

data : array_like (N, ) or (N, N)

The data representing the covariance matrix. If 1D, the input is interpreted as the diagonal elements, with all other elements zero

k_center : array_like

The wavenumbers where the center of the measurement k-bins are defined, for each element of the data matrix along one axis

mu_center : array_like,

The values where the center of the measurement mu-bins are defined, for each element of the data matrix along one axis

classmethod from_plaintext(filename)

Load the covariance matrix from a plain text file

classmethod periodic_gaussian_covariance(model, k, mu_edges, nbar, volume, Nmu=100)

Return the Gaussian prediction for the covariance between \(P(k,\mu)\) wedges for a periodic box simulation, where the number density is constant.

See eq. 17 of Grieb et al. 2015 arxiv:1509.04293

Parameters:

model : GalaxySpectrum, QuasarSpectrum

the model instance used to evaluate the theoretical \(P(k,\mu)\)

k : array_like

the array of wavenumbers (units of \(h/\mathrm{Mpc}\)) where the covariance matrix will be evaluated

mu_edges : array_like

the edges of the \(\mu\) bins

nbar : float

the constant number density in the box in units of \((\mathrm{Mpc}/h)^{-3}\) – the shot noise contribution to the covariance is the inverse of this value

volume : float

the volume of the box, in units of \((\mathrm{Mpc}/h)^3\)

Nmu : int, optional

the number of mu bins to use when performing the multipole integration over \(\mu\)

Returns:

PkmuCovarianceMatrix :

the covariance matrix object holding the Gaussian prediction for the covariance between the specified wedges

Examples

>>> import numpy
>>> from pyRSD.rsd import GalaxySpectrum
>>> from pyRSD.rsdfit.data import PkmuCovarianceMatrix
>>> volume = 1380.0**3
>>> nbar = 3e-4
>>> model = GalaxySpectrum(params='boss_dr12_fidcosmo.ini')
>>> k = numpy.arange(0., 0.4, 0.005) + 0.005/2
>>> mu_edges = [0., 0.2, 0.4, 0.6, 0.8, 1.0]
>>> C = PkmuCovarianceMatrix.periodic_gaussian_covariance(model, k, ells, nbar, volume)
to_plaintext(filename)

Save the covariance matrix to a plain text file