UFalcon package

Submodules

UFalcon.bessel module

UFalcon.bessel.get_der_jl(x, ell)[source]

Returns the spherical bessel function derivative at x for mode ell.

Parameters:
  • x (float) – x coordinate.

  • ell (int) – Desired spherical Bessel mode.

Returns:

Derivative of the Spherical Bessel Function value at x for mode ell.

Return type:

float

UFalcon.bessel.get_der_qln(ell, nmax, nstop=100, zerolminus2=None, zerolminus1=None)[source]

Returns the zeros of the spherical Bessel function derivative. Begins by assuming that the zeros of the spherical Bessel function for l lie exactly between the zeros of the Bessel function between l and l+1. This allows us to use scipy’s jn_zeros function. However, this function fails to return for high n. To work around this we estimate the first 100 zeros using scipy’s jn_zero function and then iteratively find the roots of the next zero by assuming the next zero occurs pi away from the last one. Brent’s method is then used to find a zero between pi/2 and 3pi/2 from the last zero.

Parameters:
  • ell (int) – Desired spherical Bessel mode.

  • nmax (int) – The maximum zero found for the spherical Bessel Function.

  • nstop (int) – For n <= nstop we use scipy’s jn_zeros to guess where the first nstop zeros are. These estimates are improved using Brent’s method and assuming zeros lie between -pi/2 and pi/2 from the estimates.

Returns:

Derivative of the Spherical Bessel Function value at x for mode ell.

Return type:

float

UFalcon.bessel.get_jl(x, ell)[source]

Returns the spherical bessel function at x for mode ell.

Parameters:
  • x (float) – x coordinate.

  • ell (int) – Desired spherical Bessel mode.

Returns:

Spherical Bessel Function value at x for mode ell.

Return type:

float

UFalcon.bessel.get_qln(ell, nmax, nstop=100, zerolminus2=None, zerolminus1=None)[source]

Returns the zeros of the spherical Bessel function. Begins by assuming that the zeros of the spherical Bessel function for l lie exactly between the zeros of the Bessel function between l and l+1. This allows us to use scipy’s jn_zeros function. However, this function fails to return for high n. To work around this we estimate the first 100 zeros using scipy’s jn_zero function and then iteratively find the roots of the next zero by assuming the next zero occurs pi away from the last one. Brent’s method is then used to find a zero between pi/2 and 3pi/2 from the last zero.

Parameters:
  • ell (int) – Desired spherical Bessel mode.

  • nmax (int) – The maximum zero found for the spherical Bessel Function.

  • nstop (int) – For n <= nstop we use scipy’s jn_zeros to guess where the first nstop zeros are. These estimates are improved using Brent’s method and assuming zeros lie between -pi/2 and pi/2 from the estimates.

Returns:

Derivative of the Spherical Bessel Function value at x for mode ell.

Return type:

float

UFalcon.constants module

UFalcon.construct_maps module

class UFalcon.construct_maps.construct_map_cosmogrid(maps, z_low, z_high, nside, boxsize, cosmo, n_particles, zi=0, zf=2)[source]

Bases: object

Class to compute cosmological signal maps given an input set of lightcone shells and some deatils on the simualtions

Parameters:
  • maps (ndarray) – An array of HEALPix maps that represent the lightcone shells from a simulation as raw projected particle counts

  • z_low (ndarray) – An array of floats corresponding to the lower redshift value of each lightcone shell

  • z_high (ndarray) – An array of floats corresponding to the upper redshift value of each lightcone shell

  • nside (int) – HEALPix nside of output maps (must be lower or equal to the nside of the input simualtion lightcone maps)

  • boxsize (float) – size of the box in Gigaparsec

  • cosmo (Astropy cosmology instance) – Astropy.Cosmo instance, controls the cosmology used

  • n_particles (int) – total number of parciles used in the N-Body simulation, i.e. (No. of part. per side)^3

  • zi (float) – the starting redshift at which lightcone maps should be considered (default 0)

  • zf (float) – the final redshift at which lightcone maps should be considered (default 2)

construct_IA_map(n_of_z, shift_nz=0.0, IA=None, eta=0.0, z_0=0.5, fast_mode=False)[source]

Computes IA only-map from precomputed shells in numpy format containing particle counts.

Parameters:
  • n_of_z (str or ndarray) – Either path to file containing n(z), assumed to be a text file readable with numpy with the first column containing z and the second column containing n(z), a 2D array of shape (N, 2) or a callable that is directly a redshift distribution.

  • shift_nz (float) – Shift the n(z) function by some fixed redshift (intended for easier implementation of photo z bias).

  • IA (float) – Intrinsic alignment amplitude for the NLA model.

  • eta (float) – Parameter for the redshift dependence of the NLA model.

  • z_0 (float) – Pivot parameter for the redshift dependence of the NLA model.

  • fast_mode (bool) – Instead of using quad from scipy, use a simple simpson rule. Note that this will drastically decrease the runtime of the weight calculation if your n(z) is not continuous, while reducing the accuracy and increasing the memory usage. This mode is intended for large scale runs that use the same or very similar redshift distributions. You should always benchmark the fast mode against the normal mode before you start such a run to see if the accuracy is sufficient. The possibility to tweak the accuracy parameter of the fast mode will come in future versions.

Returns:

HEALPix array representing IA map with shape (Npix,).

Return type:

ndarray

construct_delta_map(n_of_z, lin_bias)[source]

Computes clustering map from precomputed shells in numpy format containing particles counts.

Parameters:
  • n_of_z (str or ndarray) – Either path to file containing n(z), assumed to be a text file readable with numpy with the first column containing z and the second column containing n(z), a 2D array of shape (N, 2) or a callable that is directly a redshift distribution.

  • lin_bias (float) – the linear bias applied to each clustering map

Returns:

HEALPix array representing clustering map with shape (Npix,).

Return type:

ndarray

construct_isw_map(temp_path='./temp/', zmax_for_box=None, alpha=None)[source]

Computes cmb ISW from precomputed shells in numpy format containing particles counts. This function is still experimental and should be used with caution. In particular, we have found unphysical features due to box-edges in some of the maps computed using this algorithm and we caution the user to visually and statistically inspect any produced maps for unphysical features. NB this only computes the contribution to the ISW signal until z=zmax_for_box which is defined as the redshift beyond which there would be significant features due to box replication.

Parameters:
  • temp_path (str) – the temporary path where intermediate alms will be stored (recommended to put this in the local scratch space)

  • zmax_for_box (float) – a maximum user-defined redshift that the ISW contribution is computed until. If not supplied the default value of zmax_for_box=1.5*z(boxsize) i.e. the redshift at the edge of a single box multiplied by 1.5. The factor of 1.5 is motivated by experimentation with boxsizes 900 Mpc/h and 2250 Mpc/h.

  • alpha (callable) – a function that takes an input redshift (alpha(z)). This determines the factor by which the nominal minimum k mode (2pi/lbox) should be multiplied by to avoid unphysical features for redshifts ranges requiring box replications. If nothing supplied the default function is used alpha(z)=1.5 + 0.02 * np.exp((z - zmax_for_box)) which has been shown to remove unphysical features for a Cosmogrid-like set of simulations.

Returns:

HEALPix array representing CMB ISW map with shape (Npix,) in units of K.

Return type:

ndarray

construct_kappa_cmb_map()[source]

Computes cmb kappa_map from precomputed shells in numpy format containing particles counts. The redshift of recombination is fixed to 1090.0 and this code prints out the zmax up to which the CMB lensing contribution has been calculated (beyond this consider adding a gaussian synfast realization).

Returns:

HEALPix array representing CMB lensing map with shape (Npix,).

Return type:

ndarray

construct_kappa_map(n_of_z, shift_nz=0.0, IA=None, eta=0.0, z_0=0.5, fast_mode=False, fast_mode_num_points_1d=13, fast_mode_num_points_2d=512)[source]

Computes kappa_map from precomputed shells in numpy format containing particles counts User can define whether they desire a pure cosmological signal or cosmogical signal + IA contribution (nb if IA is not 0 then returns signal + IA contribution)

Parameters:
  • n_of_z (str or ndarray) – Either path to file containing n(z), assumed to be a text file readable with numpy with the first column containing z and the second column containing n(z), a 2D array of shape (N, 2) or a callable that is directly a redshift distribution.

  • shift_nz (float) – Shift the n(z) function by some fixed redshift (intended for easier implementation of photo z bias).

  • IA (float) – Intrinsic alignment amplitude for the NLA model.

  • eta (float) – Parameter for the redshift dependence of the NLA model.

  • z_0 (float) – Pivot parameter for the redshift dependence of the NLA model.

  • fast_mode – Instead of using quad from scipy, use a simple simpson rule, note that this will drastically decrease the runtime of the weight calculation if you n(z) is not continuous, while reducing the accuracy and increasing the memory usage. This mode is intended for large scale runs that use the same or very similar redshift distributions. You should always benchmark the fast mode against the normal mode before you start such a run to see if the accuracy is sufficient. The possibility to tweak the accuracy parameter of the fast mode will come in future versions.

  • fast_mode_num_points_1d (int) – number of points for computing 1D integral for weak lensing fast mode (default is tested for cosmogrid lightcones and KiDS-1000 set-up)

  • fast_mode_num_points_2d (int) – number of points for computing 2D integral for weak lensing fast mode (default is tested for cosmogrid lightcones and KiDS-1000 set-up)

Returns:

HEALPix array representing kappa map with shape (Npix,).

Return type:

ndarray

get_SBT_class_instance()[source]

Function to initialize the SBT class for the ISW map making algorithm.

get_zero_mean_overdensity_from_map(map_slice)[source]

Function to get the zero-mean overdensity given an input lightcone shell

Parameters:

map_slice (ndarray) – lightcone shell

Returns:

zero mean overdensity map

Return type:

ndarray

UFalcon.probe_weights module

class UFalcon.probe_weights.Continuous(n_of_z, interpolation_kind='linear', z_lim_low=0, z_lim_up=None, shift_nz=0.0)[source]

Bases: object

get_weight_norm(z_low, z_up, cosmo)[source]

Computes the normalization for lensing integrals.

Parameters:
  • z_low (float) – lower end of the redshift interval

  • z_up (float) – upper end of the redshift interval

class UFalcon.probe_weights.Continuous_clustering(n_of_z, interpolation_kind='linear', z_lim_low=0, z_lim_up=None, shift_nz=0.0)[source]

Bases: Continuous

Computes the clustering weights for a continuous, user-defined n(z) distribution.

class UFalcon.probe_weights.Continuous_intrinsic_alignment(n_of_z, interpolation_kind='linear', z_lim_low=0, z_lim_up=None, shift_nz=0.0, IA=1.0, eta=0.0, z_0=0.5)[source]

Bases: Continuous

Computes the intrinsic alignment weights for a continuous, user-defined n(z) distribution. N.B: the weight should then be multiplied by delta_prefactor.

class UFalcon.probe_weights.Continuous_lensing(n_of_z, interpolation_kind='linear', z_lim_low=0, z_lim_up=None, shift_nz=0.0, fast_mode=False, fast_mode_num_points_1d=13, fast_mode_num_points_2d=512)[source]

Bases: Continuous

Computes the lensing weights for a continuous, user-defined n(z) distribution. The weight should then be multiplied by kappa_prefactor.

class UFalcon.probe_weights.Dirac_lensing(z_source)[source]

Bases: object

Computes the lensing weights for a single-source redshift.

UFalcon.probe_weights.F_NLA_model(z, IA, eta, z_0, cosmo)[source]

Calculates the NLA kernel used to calculate the IA shell weight.

Parameters:
  • z (float) – Redshift where to evaluate

  • IA (float) – Galaxy intrinsic alignment amplitude

  • eta (float) – Galaxy Intrinsic alignment redshift dependence

  • z_0 (float) – Pivot parameter for the redshift dependence of the NLA model

  • cosmo (Astropy cosmology instance) – Astropy.Cosmo instance, Astropy.Cosmo instance, controls the cosmology used

Returns:

NLA kernel at redshift z

Return type:

float

UFalcon.probe_weights.delta_prefactor(n_pix, n_particles, boxsize, cosmo)[source]

Computes the prefactor to transform from number of particles to clustering see: https://arxiv.org/pdf/2007.05735.pdf.

Parameters:
  • n_pix (int) – number of healpix pixels used

  • n_particles (int) – number of particles

  • boxsize (float) – size of the box in Gigaparsec

  • cosmo (Astropy cosmology instance) – Astropy.Cosmo instance, Astropy.Cosmo instance, controls the cosmology used

Returns:

clustering prefactor

Return type:

float

UFalcon.probe_weights.delta_prefactor_nosim(cosmo)[source]

Computes the simualtion independent prefactor to transform to clustering see: https://arxiv.org/pdf/2007.05735.pdf, cosmology dependent part. This function does not include the factor depending on box size, number of particles and number of pixels. Remember to also multiply by this using the function sim_spec_prefactor.

Parameters:

cosmo (Astropy cosmology instance) – Astropy.Cosmo instance, controls the cosmology used

Returns:

cosmology dependent clustering prefactor

Return type:

float

UFalcon.probe_weights.kappa_prefactor(n_pix, n_particles, boxsize, cosmo)[source]

Computes the prefactor to transform from number of particles to convergence, see https://arxiv.org/abs/0807.3651, eq. (A.1).

Parameters:
  • n_pix (int) – number of healpix pixels used

  • n_particles (int) – number of particles

  • boxsize (float) – size of the box in Gigaparsec

  • cosmo (Astropy cosmology instance) – Astropy.Cosmo instance, Astropy.Cosmo instance, controls the cosmology used

Returns:

convergence prefactor

Return type:

float

UFalcon.probe_weights.kappa_prefactor_nosim(cosmo)[source]

Computes the simualtion independent prefactor to transform to convergence, see https://arxiv.org/abs/0807.3651, eq. (A.1), cosmology dependent part. This function does not include the factor depending on box size, number of particles and number of pixels. Remember to also multiply by this using the function sim_spec_prefactor.

Parameters:

cosmo (Astropy cosmology instance) – Astropy.Cosmo instance, controls the cosmology used

Returns:

cosmology dependent convergence prefactor

Return type:

float

UFalcon.probe_weights.sim_spec_prefactor(n_pix, box_size, n_particles)[source]

Calculate the prefactors depending on the box size, number of particles and number of pixels. See https://arxiv.org/abs/0807.3651, https://arxiv.org/pdf/2007.05735.pdf. This should be used together with delta_prefactor and kappa_prefactor.

Parameters:
  • n_pix (int) – number of healpix pixels used

  • boxsize (float) – size of the box in Gigaparsec

  • n_particles (int) – number of particles

Return prefactor:

prefactor to use

Return type:

float

UFalcon.probe_weights.w_IA(IA, eta, z_low, z_up, cosmo, nz_intpt, z_0=0.5, points=None)[source]

Calculates the weight per slice for the NLA model given a distribution of source redshifts n(z).

Parameters:
  • IA (float) – Galaxy intrinsic alignments amplitude

  • eta (float) – Galaxy Intrinsic alignment redshift dependence

  • z_low (float) – Lower redshift limit of the shell

  • z_up (float) – Upper redshift limit of the shell

  • cosmo (Astropy cosmology instance) – Astropy.Cosmo instance, Astropy.Cosmo instance, controls the cosmology used

  • nz_intpt – nz function

  • z_0 (float) – Pivot parameter for the redshift dependence of the NLA model

  • points (ndarray) – array-like, Points in redshift where integrand is evaluated (used for better numerical integration), can be None

Returns:

Shell weight for NLA model

Return type:

float

UFalcon.sbt_class module

class UFalcon.sbt_class.SphericalBesselISW(cosmo)[source]

Bases: object

Modified version of https://github.com/knaidoo29/pyGenISW/. Class for computing the ISW using spherical Bessel Transforms from maps of the density contrast given in redshift slices.

Parameters:

cosmo (Astropy cosmology instance) – Astropy.Cosmo instance, controls the cosmology used when computing linear growth

alm2sbt()[source]

Converts spherical harmonic coefficients in redshift slices to spherical Bessel coefficients. Stored as delta_lmn in units of (Mpc/h)^(1.5).

Returns:

None

Return type:

None

calc_table(zmin=0.0, zmax=10.0, zbin_num=1000, zbin_mode='linear', alpha=0.55, kind='cubic')[source]

Constructs table of cosmological linear functions to be interpolated for speed.

Parameters:
  • zmin (float) – Minimum redshift for tabulated values of the linear growth functions.

  • zmax (float) – Maximum redshift for tabulated values of the linear growth functions.

  • zbin_num (int) – Number of redshift values to compute the growth functions.

  • zbin_mode (str) – Redshift binning, either linear or log of 1+z.

  • alpha (float) – The power in the approximation to f(z) = Omega_m(z)**alpha.

  • kind (str) – The kind of interpolation used by the created interpolation functions as function of z and r.

Returns:

None

Return type:

None

create_folder(root, path=None)[source]

Creates a folder with the name ‘root’ either in the current folder if path is None or a specified path. :param root: The name of the created folder. :type root: str :param path: The name of the path of the created folder. :type path: str, optional :return: None :rtype: None

get_Dr(r)[source]

Compute the linear growth factor D at the input comoving distance r, using interpolation.

Parameters:

r (float) – Input comoving distance value.

Returns:

Corresponding linear growth factor.

Return type:

float

get_Dz(z)[source]

Compute the linear growth factor D at the input redshift z, using a direct calculation with this instance’s cosmology object.

Parameters:

z (float) – Input redshift value.

Returns:

Corresponding linear growth factor.

Return type:

float

get_Hr(r)[source]

Compute the Hubble parameter H at the input comoving distance r, using interpolation.

Parameters:

r (float) – Input comoving distance value.

Returns:

Corresponding Hubble parameter.

Return type:

float

get_Hz(z)[source]

Compute the Hubble parameter H at the input redshift z, using a direct calculation with this instance’s cosmology object.

Parameters:

z (float) – Input redshift value.

Returns:

Corresponding Hubble parameter.

Return type:

float

get_a(r)[source]

Compute the scale factor a as a function of comoving distance r from an interpolation. Wrapper around utils function.

Parameters:

r (float) – Comoving distance in Mpc

Returns:

scale factor

Return type:

float

get_fr(r)[source]

Compute the linear growth rate f at the input comoving distance r, using interpolation.

Parameters:

r (float) – Input comoving distance value.

Returns:

Corresponding linear growth rate.

Return type:

float

get_fz(z)[source]

Compute the linear growth rate f at the input redshift z, using interpolation.

Parameters:

z (float) – Input redshift value.

Returns:

Corresponding linear growth rate.

Return type:

float

get_fz_numerical(z, Dz, **kwargs)[source]

Compute the linear growth rate f at the input redshifts z, by numerically differentiating the input linear growth factors Dz.

Parameters:
  • z (array) – Input redshift value.

  • Dz (array) – Input linear growth factor value.

Returns:

Corresponding linear growth rates.

Return type:

array

get_rz(z)[source]

Compute the corresponding comoving distance value for the input redshift, using this instance’s linear interpolator.

Parameters:

z (float) – Input redshift value.

Returns:

Corresponding comoving distance value.

Return type:

float

get_rz_no_interp(z)[source]

Compute the corresponding comoving distance value for the input redshift, without using interpolation. Instead, use the cosmology instance for a direct calculation.

Parameters:

z (float) – Input redshift value.

Returns:

Corresponding comoving distance value.

Return type:

float

get_zr(r)[source]

Compute the corresponding redshift for the input comoving distance, using interpolation.

Parameters:

r (float) – Input comoving distance value.

Returns:

Corresponding redshift value.

Return type:

float

get_zr_no_interp(r)[source]

Compute the corresponding redshift for the input comoving distance, without using interpolation. Instead, use the cosmology instance for a direct calculation.

Parameters:

r (float) – Input comoving distance value.

Returns:

Corresponding redshift value.

Return type:

float

numerical_differentiate(x, f, equal_spacing=False, interpgrid=1000, kind='cubic')[source]

For unequally spaced data we interpolate onto an equal spaced 1d grid. Then, we use the symmetric two-point derivative and the non-symmetric three point derivative estimator.

Parameters:
  • x (array) – X-axis.

  • f (array) – Function values at x.

  • equal_spacing (bool, optional) – Automatically assumes data is not equally spaced and will interpolate from it.

  • interpgrid (int, optional) – Grid spacing for the interpolation grid, if equal spacing is False.

  • kind (str, optional) – Interpolation kind.

Returns:

df, the numerical differentiation values for f evaluated at points x.

Return type:

array

Notes

For non-boundary values: df/dx = (f(x + dx) - f(x - dx)) / (2*dx)

For boundary values: df/dx = (- f(x + 2dx) + 4f(x + dx) - 3f(x)) / (2*dx)

save_sbt(prefix=None)[source]

Saves spherical Bessel transform coefficients.

Parameters:

prefix (str) – Prefix for file containing spherical Bessel transform.

Returns:

None

Return type:

None

sbt2isw_alm(zmin=None, zmax=None)[source]

Compute the ISW spherical harmonics between zmin and zmax from the computed spherical Bessel Transform.

Parameters:
  • zmin (float) – Minimum redshift for ISW computation.

  • zmax (float) – Maximum redshift for ISW computation.

Returns:

Array of alm coefficients for the spherical harmonics.

Return type:

array

sbt2isw_map(zmin, zmax, nside=256)[source]

Returns a healpix map of the ISW between zmin and zmax computed from the spherical Bessel Transform.

Parameters:
  • zmin (float) – Minimum redshift for ISW computation.

  • zmax (float) – Maximum redshift for ISW computation.

  • nside (int) – NSIDE for healpix map.

Returns:

Corresponding ISW map.

Return type:

array

setup(zmin, zmax, zedge_min, zedge_max, temp_path, alpha, kmin=None, kmax=0.1, lmax=None, nmax=None, uselightcone=True, boundary_conditions='normal')[source]

Configure instance for future computations.

Parameters:
  • zmin (float) – Minimum redshift for spherical Bessel transform.

  • zmax (float) – Maximum redshift for spherical Bessel transform.

  • zedge_min (array) – Minimum redshift edge for each slice.

  • zedge_max (array) – Maximum redshift edge for each slice.

  • kmin (float) – Minium Fourier mode to consider.

  • kmax (float) – Maximum Fourier mode to consider.

  • lmax (int) – Maximum l mode to compute to, if None will be computed based on kmax.

  • nmax (int) – Maximum n mode to compute to, if None will be computed based on kmax.

  • uselightcone (bool) – True if density contrast maps are given as a lightcone and not all at redshift 0.

  • boundary_conditions (str) – Setting for boundary conditions used. This can assume the following values: - normal: boundaries where spherical bessel function is zero. - derivative: boundaries where the derivative of the spherical Bessel function is zero.

Returns:

None

Return type:

None

slice2alm(map_slice, index)[source]

Given a density contrast map and its corresponding index (used to obtain the map’s zedges minimum and maximum), slice2alm will convert the map to its spherical harmonics and save the corresponding output files.

Parameters:
  • map_slice (array) – Healpix density contrast map.

  • index (int) – Index of the slice, used to determine the slice’s zedges.

Returns:

None

Return type:

None

z2a(z)[source]

Convert an input redshift z into the corresponding scale factor a.

Parameters:

z (float or array) – Input redshift value.

Returns:

Corresponding scale factor.

Return type:

float or array

UFalcon.shells module

UFalcon.shells.construct_shells(dirpath, z_shells, boxsize, cosmo, nside, file_format='l-picola')[source]

Reads in particle positions stored in all the binary file produced by either L-PICOLA or PKDGRAV and transforms their angular positions to counts in healpix pixels corresponding to shells at different redshifts.

Parameters:
  • dirpath (str) – path to the directory holding the binary files with particle positions

  • z_shells (ndarray) – array containing the discrete redshifts steps, over which the lightcone is constructed

  • boxsize (float) – size of the box in Gigaparsec

  • cosmo (Astropy cosmology instance) – Astropy.Cosmo instance, controls the cosmology used

  • nside (int) – nside of the healpix map

  • file_format (str) – data format, either l-picola or pkdgrav

Returns:

matrix with dimension (len(z_shells) - 1, Npix) containing the number counts for each shell and pixel-index

Return type:

ndarray

UFalcon.shells.read_file(path, boxsize, cosmo, file_format='l-picola')[source]

Reads in particle positions stored in a binary file produced by either L-PICOLA or PKDGRAV.

Parameters:
  • path (str) – path to binary file holding particle positions

  • boxsize (float) – size of the box in Gigaparsec

  • cosmo (Astropy cosmology instance) – Astropy.Cosmo instance, controls the cosmology used

  • file_format (str) – data format, either l-picola or pkdgrav

Returns:

theta- and phi-coordinates of particles inside the shell

Return type:

tuple

UFalcon.shells.read_lpicola(path, h, boxsize)[source]

Reads in a binary data file produced by L-Picola.

Parameters:
  • path (str) – path to file

  • h (float) – Dimensionless Hubble parameter

  • boxsize (float) – size of the box in Gigaparsec

Returns:

3-tuple containing (x, y, z) particle positions in Megaparsec

Return type:

tuple

UFalcon.shells.read_pkdgrav(path, boxsize, n_rows_per_block=1000000)[source]

Reads in a binary data file produced by PKDGRAV.

Parameters:
  • path (str) – path to file

  • boxsize (float) – size of the box in Gigaparsec

  • n_rows_per_block (int) – number of rows to read in one block, allows to limit memory consumption for large files

Returns:

3-tuple containing (x, y, z) particle positions in Megaparsec

Return type:

tuple

UFalcon.shells.thetaphi_to_pixelcounts(theta, phi, nside)[source]

Transforms angular particle positions to counts in healpix pixels. The size of the output array equals the index of the last non-empty pixel (i.e. the largest healpix index with at least one count).

Parameters:
  • theta (float) – healpix theta-coordinate

  • phi (float) – healpix phi-coordinate

  • nside (int) – nside of the healpix map

Returns:

counts in each pixel, maximum size: nside - 1

Return type:

ndarray

UFalcon.shells.xyz_to_spherical(xyz_coord)[source]

Transform from comoving cartesian (x, y, z)- to spherical coordinates (comoving radius, healpix theta, healpix phi).

Parameters:

xyz_coord – cartesian coordinates, shape: (number of particles, 3)

Returns:

comoving radius, theta, phi

Return type:

ndarray

UFalcon.utils module

UFalcon.utils.a_of_r(r, cosmo)[source]

Computes the scale factor at a given comoving distance in MpC. Scalar or vector input. Note only supports calculation up to redshift z=12.

Parameters:
  • r (float) – input comoving distance in MpC at which scale factor should we found

  • cosmo (Astropy cosmology instance) – Astropy.Cosmo instance, controls the cosmology used

Returns:

scale factor

Return type:

float

UFalcon.utils.comoving_distance(z_low, z_up, cosmo)[source]

Computes the comoving distance between two redshifts. Scalar input only.

Parameters:
  • z_low (ndarray) – lower redshift

  • z_up (ndarray) – upper redshift, must have same shape as z_low

  • cosmo (Astropy cosmology instance) – Astropy.Cosmo instance, controls the cosmology used

Returns:

comoving distance

Return type:

ndarray

UFalcon.utils.dimensionless_comoving_distance(z_low, z_up, cosmo)[source]

Computes the dimensionless comoving distance between two redshifts. Scalar input only. Legacy code - see updated function

Parameters:
  • z_low (ndarray) – lower redshift

  • z_up (ndarray) – upper redshift, must have same shape as z_low

  • cosmo (Astropy cosmology instance) – Astropy.Cosmo instance, controls the cosmology used

Returns:

dimensionless comoving distance

Return type:

ndarray

UFalcon.utils.dimensionless_comoving_distance_old(z_low, z_up, cosmo, fast_mode=0, z_max_interp=10)[source]

Computes the dimensionless comoving distance between two redshifts. Scalar input only. Legacy code - see updated function

Parameters:
  • z_low (ndarray) – lower redshift

  • z_up (ndarray) – upper redshift, must have same shape as z_low

  • cosmo (Astropy cosmology instance) – Astropy.Cosmo instance, controls the cosmology used

  • fast_mode (bool) – Instead of using quad from scipy, use a simple romberg integration rule (this works here because we know that the dimless com behaves and is differentiable)

Returns:

dimensionless comoving distance

Return type:

ndarray

UFalcon.utils.growth_z(z, cosmo)[source]

Computes the normalized linear growth factor as a function of redshift (i.e. D(z)). Scalar input only. Normalized means D(z)=1 at z=0.

Parameters:
  • z (ndarray) – redshift or array of redshifts at which growth factor should be calculated

  • cosmo (Astropy cosmology instance) – Astropy.Cosmo instance, controls the cosmology used

Returns:

linear growth factor at z i.e. D(z)

Return type:

ndarray

UFalcon.utils.kappa_to_gamma(kappa_map, lmax=None, use_pixel_weights=True)[source]

Computes a gamma_1- and gamma_2-map from a kappa-map, s.t. the kappa TT-spectrum equals the gamma EE-spectrum.

Parameters:
  • kappa_map (ndarray) – kappa map as a HEALPix array

  • lmax (float) – maximum multipole to consider, default: 3 * nside - 1

  • use_pixel_weights (bool) – Use pixelweights for the map2alm transformation. This delivers the most accurate transform according to healpy, but requires the pixel weights, which will be downloaded automatically.

Returns:

gamma_1- and gamma_2-map

Return type:

ndarray

Module contents