Various analysis functions

Derived quantities and basic statistics

c2raytools.temperature.calc_dt(xfrac, dens, z=-1)

Calculate the differential brightness temperature assuming T_s >> T_CMB

Parameters:
  • xfrac (XfracFile object, string or numpy array): the ionization fraction
  • dens (DensityFile object, string or numpy array): density in cgs units
  • z = -1 (float): The redshift (if < 0 this will be figured out from the files)
Returns:
The differential brightness temperature as a numpy array with the same dimensions as xfrac.
c2raytools.temperature.calc_dt_lightcone(xfrac, dens, lowest_z, los_axis=2)

Calculate the differential brightness temperature assuming T_s >> T_CMB for lightcone data.

Parameters:
  • xfrac (string or numpy array): the name of the ionization

    fraction file (must be cbin), or the xfrac lightcone data

  • dens (string or numpy array): the name of the density

    file (must be cbin), or the density data

  • lowest_z (float): the lowest redshift of the lightcone volume

  • los_axis = 2 (int): the line-of-sight axis

Returns:
The differential brightness temperature as a numpy array with the same dimensions as xfrac.
c2raytools.temperature.mean_dt(z)

Get the mean dT at redshift z

Parameters:
  • z (float or numpy array): the redshift
Returns:
dT (float or numpy array) the mean brightness temperature in mK
c2raytools.statistics.apply_func_along_los(signal, func, los_axis)

Apply a function, such as np.var() or np.mean(), along the line-of-sight axis of a signal on a per-slice basis.

Parameters:
  • signal (numpy array): the signal
  • func (callable): the function to apply
  • los_axis (integer): the line-of-sight axis
Returns:
An array of length signal.shape[los_axis]
Example:

Calculate the variance of a lightcone along the line-of-sight:

>>> lightcone = c2t.read_cbin('my_lightcone.cbin')
>>> dT_var = c2t.apply_func_along_los(lightcone, np.var, 2)
c2raytools.statistics.mass_weighted_mean_xi(xi, rho)

Calculate the mass-weighted mean ionization fraction.

Parameters:
  • xi (numpy array): the ionized fraction
  • rho (numpy array): the density (arbitrary units)
Returns:
The mean mass-weighted ionized fraction.
c2raytools.statistics.signal_overdensity(signal, los_axis)

Divide by the mean of the signal along the los axis and subtract one.

Parameters:
  • signal (numpy array): the signal to subtract the mean from
  • los_axis (integer): the line-of-sight axis
Returns:
The signal with the mean subtracted

TODO:vectorize

c2raytools.statistics.skewness(x)

Calculate the skewness of an array. Note that IDL calculates the skewness in a slightly different way than Python. This routine uses the IDL definition.

Parameters:
  • x (numpy array): The array containing the input data
Returns:
The skewness.
c2raytools.statistics.subtract_mean_signal(signal, los_axis)

Subtract the mean of the signal along the los axis.

Parameters:
  • signal (numpy array): the signal to subtract the mean from
  • los_axis (integer): the line-of-sight axis
Returns:
The signal with the mean subtracted

TODO:vectorize

c2raytools.tau.tau(ionfractions, redshifts, num_points=50)

Calculate the optical depth to Thomson scattering.

Parameters:
  • ionfractions (numpy array): an array containing the ionized fraction

    in various points along the line-of-sight.

  • redshifts (numpy array): an array containing the redshifts of

    the same points as ionfractions. Must be the same length as ionfractions

  • num_points = 50 (integer): the number of points used for the integration

Returns:

Tuple containing (tau_0, tau_z)

tau_0 is the optical depth at each redshift

tau_z is the corresponding redshift

Notes:
  • The Universe is assumed to be fully ionized at the lowest redshift supplied.
  • To get the total optical depth, look at the last value in tau_0
Example:

To calculate the optical depth for a scenario where the Universe is instantaneously reionized:

>>> z_reion = 11.
>>> redshifts = np.linspace(z_reion, 1100., 50)
>>> ionfracs = np.zeros(len(redshifts))
>>> tau0, tau_z = tau(ionfracs, redshifts)
>>> print 'Total tau: ', tau0[-1]
0.0884755058758

Power spectra

c2raytools.power_spectrum.cross_power_spectrum_1d(input_array1_nd, input_array2_nd, kbins=100, box_dims=None, return_n_modes=False)

Calculate the spherically averaged cross power spectrum of two arrays and return it as a one-dimensional array.

Parameters:
  • input_array1_nd (numpy array): the first data array

  • input_array2_nd (numpy array): the second data array

  • kbins = 100 (integer or array-like): The number of bins,

    or a list containing the bin edges. If an integer is given, the bins are logarithmically spaced.

  • box_dims = None (float or array-like): the dimensions of the

    box. If this is None, the current box volume is used along all dimensions. If it is a float, this is taken as the box length along all dimensions. If it is an array-like, the elements are taken as the box length along each axis.

  • return_n_modes = False (bool): if true, also return the

    number of modes in each bin

Returns:
A tuple with (Pk, bins), where Pk is an array with the cross power spectrum and bins is an array with the k bin centers.
c2raytools.power_spectrum.cross_power_spectrum_mu(input_array1, input_array2, los_axis=0, mubins=20, kbins=10, box_dims=None, weights=None, exclude_zero_modes=True)

Calculate the cross power spectrum and bin it in mu=cos(theta) and k input_array is the array to calculate the power spectrum from

Parameters:
  • input_array1 (numpy array): the first data array

  • input_array2 (numpy array): the second data array

  • los_axis = 0 (integer): the line-of-sight axis

  • mubins = 20 (integer): the number of mu bins

  • kbins = 10 (integer or array-like): The number of bins,

    or a list containing the bin edges. If an integer is given, the bins are logarithmically spaced.

  • box_dims = None (float or array-like): the dimensions of the

    box. If this is None, the current box volume is used along all dimensions. If it is a float, this is taken as the box length along all dimensions. If it is an array-like, the elements are taken as the box length along each axis.

  • exlude_zero_modes = True (bool): if true, modes with any components

    of k equal to zero will be excluded.

Returns:
A tuple with (Pk, mubins, kbins), where Pk is an array with the cross power spectrum of dimensions (n_mubins x n_kbins), mubins is an array with the mu bin centers and kbins is an array with the k bin centers.
TODO:
Add support for (non-numpy) lists for the bins
c2raytools.power_spectrum.cross_power_spectrum_nd(input_array1, input_array2, box_dims)

Calculate the cross power spectrum two arrays and return it as an n-dimensional array, where n is the number of dimensions in input_array box_side is the size of the box in comoving Mpc. If this is set to None (default), the internal box size is used

Parameters:
  • input_array1 (numpy array): the first array to calculate the

    power spectrum of. Can be of any dimensions.

  • input_array2 (numpy array): the second array. Must have same

    dimensions as input_array1.

  • box_dims = None (float or array-like): the dimensions of the

    box. If this is None, the current box volume is used along all dimensions. If it is a float, this is taken as the box length along all dimensions. If it is an array-like, the elements are taken as the box length along each axis.

Returns:
The cross power spectrum in the same dimensions as the input arrays.
TODO:
Also return k values.
c2raytools.power_spectrum.mu_binning(powerspectrum, los_axis=0, mubins=20, kbins=10, box_dims=None, weights=None, exclude_zero_modes=True)

This function is for internal use only.

c2raytools.power_spectrum.power_spectrum_1d(input_array_nd, kbins=100, box_dims=None, return_n_modes=False)

Calculate the spherically averaged power spectrum of an array and return it as a one-dimensional array.

Parameters:
  • input_array_nd (numpy array): the data array

  • kbins = 100 (integer or array-like): The number of bins,

    or a list containing the bin edges. If an integer is given, the bins are logarithmically spaced.

  • box_dims = None (float or array-like): the dimensions of the

    box. If this is None, the current box volume is used along all dimensions. If it is a float, this is taken as the box length along all dimensions. If it is an array-like, the elements are taken as the box length along each axis.

  • return_n_modes = False (bool): if true, also return the

    number of modes in each bin

Returns:
A tuple with (Pk, bins), where Pk is an array with the power spectrum and bins is an array with the k bin centers.
c2raytools.power_spectrum.power_spectrum_mu(input_array, los_axis=0, mubins=20, kbins=10, box_dims=None, weights=None, exclude_zero_modes=True)

Calculate the power spectrum and bin it in mu=cos(theta) and k input_array is the array to calculate the power spectrum from

Parameters:
  • input_array (numpy array): the data array

  • los_axis = 0 (integer): the line-of-sight axis

  • mubins = 20 (integer): the number of mu bins

  • kbins = 10 (integer or array-like): The number of bins,

    or a list containing the bin edges. If an integer is given, the bins are logarithmically spaced.

  • box_dims = None (float or array-like): the dimensions of the

    box. If this is None, the current box volume is used along all dimensions. If it is a float, this is taken as the box length along all dimensions. If it is an array-like, the elements are taken as the box length along each axis.

  • exlude_zero_modes = True (bool): if true, modes with any components

    of k equal to zero will be excluded.

Returns:
A tuple with (Pk, mubins, kbins), where Pk is an array with the power spectrum of dimensions (n_mubins x n_kbins), mubins is an array with the mu bin centers and kbins is an array with the k bin centers.
c2raytools.power_spectrum.power_spectrum_nd(input_array, box_dims=None)

Calculate the power spectrum of input_array and return it as an n-dimensional array, where n is the number of dimensions in input_array box_side is the size of the box in comoving Mpc. If this is set to None (default), the internal box size is used

Parameters:
  • input_array (numpy array): the array to calculate the

    power spectrum of. Can be of any dimensions.

  • box_dims = None (float or array-like): the dimensions of the

    box. If this is None, the current box volume is used along all dimensions. If it is a float, this is taken as the box length along all dimensions. If it is an array-like, the elements are taken as the box length along each axis.

Returns:
The power spectrum in the same dimensions as the input array.
c2raytools.power_spectrum.radial_average(input_array, box_dims, kbins=10)

Radially average data. Mostly for internal use.

Parameters:
  • input_array (numpy array): the data array

  • box_dims = None (float or array-like): the dimensions of the

    box. If this is None, the current box volume is used along all dimensions. If it is a float, this is taken as the box length along all dimensions. If it is an array-like, the elements are taken as the box length along each axis.

  • kbins = 10 (integer or array-like): The number of bins,

    or a list containing the bin edges. If an integer is given, the bins are logarithmically spaced.

Returns:
A tuple with (data, bins, n_modes), where data is an array with the averaged data, bins is an array with the bin centers and n_modes is the number of modes in each bin
c2raytools.power_legendre.power_spectrum_multipoles(input_array, kbins=10, box_dims=None, los_axis=0, output=['P0', 'P2', 'P4', 'nmodes'], exclude_zero_modes=False)

Calculate the power spectrum of an array and expand it in the first three Legendre polynomials.

Parameters:
  • input_array (numpy array): the array to calculate the

    power spectrum of. Can be of any dimensions.

  • kbins = 10 (integer or array-like): The number of bins,

    or a list containing the bin edges. If an integer is given, the bins are logarithmically spaced.

  • box_dims = None (float or array-like): the dimensions of the

    box. If this is None, the current box volume is used along all dimensions. If it is a float, this is taken as the box length along all dimensions. If it is an array-like, the elements are taken as the box length along each axis.

  • los_axis = 0 (integer): the line-of-sight axis

  • output = [‘P0’, ‘P2’, ‘P4’, ‘nmodes’] (list): the multipole moments to

    include in the output. For example, to get only the P2 moment, pass in [‘P2’]. Can also contain ‘nmodes’ to return the number of Fourier modes per bin

  • exlude_zero_modes = True (bool): if true, modes with any components

    of k equal to zero will be excluded.

Returns:
A tuple with (multipoles, k) where multipoles is a dictionary containing the multipoles (the keys are the values passed to the output parameter) and k contains the midpoints of the k bins. All arrays have the same dimension

Line-of-sight effects

c2raytools.pv_mpm.get_distorted_dt(dT, kms, redsh, los_axis=0, velocity_axis=0, num_particles=10, periodic=True)

Apply peculiar velocity distortions to a differential temperature box, using the Mesh-Particle-Mesh method, as described in http://arxiv.org/abs/1303.5627

Parameters:
  • dT (numpy array): the differential temperature box

  • kms (numpy array): velocity in km/s, array of dimensions

    (3,mx,my,mz) where (mx,my,mz) is dimensions of dT

  • redsh (float): the redshift

  • los_axis = 0 (int): the line-of-sight axis of the output volume

    (must be 0, 1 or 2)

  • velocity_axis = 0 (int): the index that indicates los velocity

  • num_particles = 10 (int): the number of particles to use per cell

    A higher number gives better accuracy, but worse performance.

  • periodic = True (bool): whether or not to apply periodic boundary

    conditions along the line-of-sight. If you are making a lightcone volume, this should be False.

Returns:
The redshift space box as a numpy array with same dimensions as dT.
Example:

Read a density file, a velocity file and an xfrac file, calculate the brightness temperature, and convert it to redshift space.

>>> vfile = c2t.VelocityFile('/path/to/data/8.515v_all.dat')
>>> dfile = c2t.DensityFile('/path/to/data/8.515n_all.dat')
>>> xfile = c2t.XfracFile('/path/to/data/xfrac3d_8.515.bin')
>>> dT = c2t.calc_dt(xfile, dfile)
>>> kms = vfile.get_kms_from_density(dfile)
>>> dT_zspace = get_distorted_dt(dT, kms, dfile.z, los_axis = 0)

Note

At the moment, it is a requirement that dimensions perpendicular to the line-of-sight are equal. For example, if the box dimensions are (mx, my, mz) and the line-of-sight is along the z axis, then mx has to be equal to my.

Note

If dT is a lightcone volume, los_axis is not necessarily the same as velocity_axis. The lightcone volume methods in c2raytools all give output volumes that have the line-of-sight as the last index, regardless of the line-of-sight axis. For these volumes, you should always use los_axis=2 and set velocity_axis equal to whatever was used when producing the real-space lightcones.

c2raytools.pv_mpm.make_pv_box(dT_filename, vel_filename, dens_filename, z, los=0, num_particles=10)

Convenience method to read files and make a distorted box.

Parameters:
  • dT_filename (string): the name of the dT file

  • vel_filename (string): the name of the velocity file

  • dens_filename (string): the name of the density file

  • z (float): the redshift

  • los (integer): the line-of-sight axis

  • num_particles (integer): the number of particles to pass

    to get_distorted_dt

Returns:
The redshift space box
c2raytools.lightcone.get_lightcone_subvolume(lightcone, redshifts, central_z, depth_mhz=None, depth_mpc=None, odd_num_cells=True, subtract_mean=True, fov_Mpc=None)

Extract a subvolume from a lightcone, at a given central redshift, and with a given depth. The depth can be specified in Mpc or MHz. You must give exactly one of these parameters.

Parameters:
  • ligthcone (numpy array): the lightcone

  • redshifts (numpy array): the redshifts along the LOS

  • central_z (float): the central redshift of the subvolume

  • depth_mhz (float): the depth of the subvolume in MHz

  • depth_mpc (float): the depth of the subvolume in Mpc

  • odd_num_cells (bool): if true, the depth of the box will always

    be an odd number of cells. This avoids problems with power spectrum calculations.

  • subtract_mean (bool): if true, subtract the mean of the signal

  • fov_Mpc (float): the FoV size in Mpc

Returns:
Tuple with (subvolume, dims) where dims is a tuple with the dimensions of the subvolume in Mpc
c2raytools.lightcone.make_lightcone(filenames, z_low=None, z_high=None, file_redshifts=None, cbin_bits=32, cbin_order='c', los_axis=0, raw_density=False, interpolation='linear')

Make a lightcone from xfrac, density or dT data. Replaces freq_box.

Parameters:
  • filenames (string or array): The coeval cubes.

    Can be either any of the following:

    • An array with the file names
    • A text file containing the file names
    • The directory containing the files (must only contain

    one type of files)

  • z_low (float): the lowest redshift. If not given, the redshift of the

    lowest-z coeval cube is used.

  • z_high (float): the highest redshift. If not given, the redshift of the

    highest-z coeval cube is used.

  • file_redshifts (string or array): The redshifts of the coeval cubes.

    Can be any of the following types:

    • None: determine the redshifts from file names
    • array: array containing the redshift of each coeval cube
    • filename: the name of a data file to read the redshifts from
  • cbin_bits (int): If the data files are in cbin format, you may specify

    the number of bits.

  • cbin_order (char): If the data files are in cbin format, you may specify

    the order of the data.

  • los_axis (int): the axis to use as line-of-sight for the coeval cubes

  • raw_density (bool): if this is true, and the data is a

    density file, the raw (simulation units) density will be returned instead of the density in cgs units

  • interpolation (string): can be ‘linear’, ‘step’, ‘sigmoid’ or

    ‘step_cell’. Determines how slices in between output redshifts are interpolated.

Returns:

(lightcone, z) tuple

lightcone is the lightcone volume where the first two axes have the same size as the input cubes

z is an array containing the redshifts along the line-of-sight

Note

If z_low is given, that redshift will be the lowest one included, even if there is no coeval box at exactly that redshift. This can give results that are subtly different from results calculated with the old freq_box routine.

c2raytools.lightcone.make_velocity_lightcone(vel_filenames, dens_filenames, z_low=None, z_high=None, file_redshifts=None, los_axis=0)

Make a lightcone from velocity data. Since velocity files contain momentum rather than actual velocity, you must specify filenames for both velocity and density.

Parameters:
  • vel_filenames (string or array): The coeval velocity cubes.

    Can be any of the following:

    • An array with the file names
    • A text file containing the file names
    • The directory containing the files (must only contain

    one type of files)

  • dens_filenames (string or array): The coeval density cubes.

    Same format as vel_filenames.

  • z_low (float): the lowest redshift. If not given, the redshift of the

    lowest-z coeval cube is used.

  • z_high (float): the highest redshift. If not given, the redshift of the

    highest-z coeval cube is used.

  • file_redshifts (string or array): The redshifts of the coeval cubes.

    Can be any of the following types:

    • None: determine the redshifts from file names
    • array: array containing the redshift of each coeval cube
    • filename: the name of a data file to read the redshifts from
  • los_axis (int): the axis to use as line-of-sight for the coeval cubes

Returns:

(lightcone, z) tuple

lightcone is the lightcone volume where the first two axes have the same size as the input cubes

z is an array containing the redshifts along the line-of-sight

c2raytools.lightcone.redshifts_at_equal_comoving_distance(z_low, z_high, box_grid_n=256, box_length_mpc=None)

Make a frequency axis vector with equal spacing in co-moving LOS coordinates. The comoving distance between each frequency will be the same as the cell size of the box.

Parameters:
  • z_low (float): The lower redshift
  • z_high (float): The upper redhisft
  • box_grid_n = 256 (int): the number of slices in an input box
  • box_length_mpc (float): the size of the box in cMpc. If None,

set to conv.LB

Returns:
numpy array containing the redshifts

Created on Sep 17, 2014

@author: Hannes Jensen

Methods to convert data between physical (cMpc) coordinates and observational (angular-frequency) coordinates

c2raytools.angular_coordinates.angular_slice_to_physical(input_slice, z, slice_size_deg, output_cell_size, output_size_mpc, order=0, prefilter=True)

Interpolate a slice in angular coordinates to physical

Returns:
(physical_slice, size_mpc)
c2raytools.angular_coordinates.bin_lightcone_in_frequency(lightcone, z_low, box_size_mpc, dnu)

Bin a lightcone in frequency bins.

Parameters:
  • lightcone (numpy array): the lightcone in length units
  • z_low (float): the lowest redshift of the lightcone
  • box_size_mpc (float): the side of the lightcone in Mpc
  • dnu (float): the width of the frequency bins in MHz
Returns:
The lightcone, binned in frequencies with high frequencies first The frequencies along the line of sight in MHz
c2raytools.angular_coordinates.bin_lightcone_in_mpc(lightcone, frequencies, cell_size_mpc)

Bin a lightcone in Mpc slices along the LoS

c2raytools.angular_coordinates.observational_lightcone_to_physical(observational_lightcone, input_freqs, input_dtheta)

Interpolate a lightcone volume measured in observational (angle/frequency) units into physical (length) units. The output resolution will be set to the coarest one, as determined either by the angular or the frequency resolution. The lightcone must have the LoS as the last index, with frequencies decreasing along the LoS.

Parameters:
  • observational_lightcone (numpy array): the input lightcone volume

  • input_freqs (numpy array): the frequency in MHz of each slice along the

    line of sight of the input

  • input_dheta (float): the angular size of a cell in arcmin

Returns:
  • The output volume
  • The redshifts along the LoS of the output
  • The output cell size in Mpc
c2raytools.angular_coordinates.physical_lightcone_to_observational(physical_lightcone, input_z_low, output_dnu, output_dtheta, input_box_size_mpc=None)

Interpolate a lightcone volume from physical (length) units to observational (angle/frequency) units.

Parameters:
  • physical_lightcone (numpy array): the lightcone volume

  • input_z_low (float): the lowest redshift of the input lightcone

  • output_dnu (float): the frequency resolution of the output volume in MHz

  • output_dtheta (float): the angular resolution of the output in arcmin

  • input_box_size_mpc (float): the size of the input FoV in Mpc.

    If None (default), this will be set to conv.LB

Returns:
  • The output volume as a numpy array
  • The output frequencies in MHz as an array of floats
c2raytools.angular_coordinates.physical_slice_to_angular(input_slice, z, slice_size_mpc, fov_deg, dtheta, order=0)

Interpolate a slice in physical coordinates to angular coordinates.

Parameters:
  • input_slice (numpy array): the 2D slice in physical coordinates

  • z (float): the redshift of the input slice

  • slice_size_Mpc (float): the size of the input slice in cMpc

  • fov_deg (float): the field-of-view in degrees. The output will be

    padded to match this size

  • dtheta (float): the target resolution in arcmin

Returns:
(angular_slice, size_deg)
c2raytools.angular_coordinates.resample_slice(input_slice, n_output_cells, order=0, prefilter=True)

Resample a 2D slice to new dimensions.

Cosmology and various conversions

c2raytools.cosmology.angular_size(dl, z)

Calculate the angular size of an object at a given redshift.

Parameters:
  • dl (float or array): the physical size in kpc
  • z (float or array): the redshift of the object
Returns:
The angluar size in arcseconds
c2raytools.cosmology.angular_size_comoving(cMpc, z)

Calculate the angular size in degrees of an object with a given comoving size.

Parameters:
  • cMpc (float or array): the size in comoving Mpc
  • z (float or array): the redshift of the object
Returns:
The angular size in degrees
c2raytools.cosmology.c_to_p(z_to_cdist, z)

Convert comoving distance to proper distance

Parameters:
  • z_to_cdist: The comoving distance
  • z: the redshift
Returns:
Proper distance
c2raytools.cosmology.cdist_to_z(dist)

Calculate the redshift correspoding to the given redshift.

Parameters:
  • dist (float or array): the distance in comoving Mpc
Returns:

redshift corresponding to the distance.

Note

Uses a precalculated table for interpolation. Only valid for 0 < z < 100

c2raytools.cosmology.deg_to_cdist(deg, z)

Calculate the size in cMpc of an object with given angular diameter.

Parameters:
  • deg (float or array): the size in degrees
  • z (float or array): the redshift
Returns:
The size in cMpc
c2raytools.cosmology.luminosity_distance(z, k=0)

Calculate the luminosity distance for a given redshift.

Parameters:
  • z (float or array): the redshift(s)
  • k = 0 (float): the curvature constant.
Returns:
The luminosity distance in Mpc
c2raytools.cosmology.nu_to_cdist(nu21)

Calculate the comoving distance to a given 21 cm frequency

Parameters:
  • nu21 (float or array): redshifted 21 cm frequency in MHz
Returns:
Comoving distance in Mpc
c2raytools.cosmology.nu_to_wavel(nu)

Convert frequency to wavelength

Parameters:
  • nu (float or array): the frequency in MHz
Returns:
The wavelength in meters
c2raytools.cosmology.nu_to_z(nu21)

Convert 21 cm frequency in MHz to redshift

Parameters:
  • nu21 (float or array): redshifted 21 cm frequency in MHz
Returns:
Redshift
c2raytools.cosmology.p_to_c(pdist, z)

Convert proper distance to comoving distance

Parameters:
  • pdist: The proper distance
  • z: the redshift
Returns:
Comoving distance
c2raytools.cosmology.z_to_cdist(z)

Calculate the comoving distance

Parameters:
z (float or array): redshift
Returns:
Comoving distance in Mpc
c2raytools.cosmology.z_to_nu(z)

Get the 21 cm frequency that corresponds to redshift z

Parameters:
  • z (float or array): redshift
Returns:
redshifted 21 cm frequency in MHz

Miscellaneous

c2raytools.beam_convolve.beam_convolve(input_array, z, fov_mpc, beam_w=None, max_baseline=None, beamshape='gaussian')

Convolve input_array with a beam of the specified form. The beam can be specified either by a width in arcminutes, or as a maximum baseline. You must specify exctly one of these parameters.

Parameters:
  • input_array (numpy array): the array to be convolved

  • z (float): the redshift of the map

  • fov_mpc (float): the field of view in cMpc

  • beam_w (float) = None: the width of the beam in arcminutes

  • max_baseline (float): the maximum baseline in meters

    (can be specified instead of beam_w)

  • beamshape (string): The shape of the beam

    (only ‘gaussian’ supported at this time)

Returns:
The convolved array (a numpy array with the same dimensions as input_array).
c2raytools.smoothing.gauss_kernel(size, sigma=1.0, fwhm=None)

Generate a normalized gaussian kernel, defined as exp(-(x^2 + y^2)/(2sigma^2)).

Parameters:
  • size (int): Width of output array in pixels.

  • sigma = 1.0 (float): The sigma parameter for the Gaussian.

  • fwhm = None (float or None): The full width at half maximum.

    If this parameter is given, it overrides sigma.

Returns:
numpy array with the Gaussian. The dimensions will be size x size or size x sizey depending on whether sizey is set. The Gaussian is normalized so that its integral is 1.
c2raytools.smoothing.get_beam_w(baseline, z)

Calculate the width of the beam for an interferometer with a given maximum baseline. It is assumed that observations are done at lambda = 21*(1+z) cm

Parameters:
  • baseline (float): the maximum baseline in meters
  • z (float): the redshift
Returns:
The beam width in arcminutes
c2raytools.smoothing.interpolate2d(input_array, x, y, order=0)

Same as interpolate2d but for 2D data

Parameters:
  • input_array (numpy array): the array to interpolate

  • x (numpy array): the output coordinates along the x axis

    expressed as (fractional) indices

  • y (numpy array): the output coordinates along the y axis

    expressed as (fractional) indices

  • order (int): the order of the spline interpolation. Default

    is 0 (linear interpolation). Setting order=1 gives the same results as IDL’s interpolate function

Returns:
Interpolated array with shape (nx, ny), where nx and ny are the lengths of the arrays x and y respectively.
c2raytools.smoothing.interpolate3d(input_array, x, y, z, order=0)

This function is a recreation of IDL’s interpolate routine. It takes an input array, and interpolates it to a new size, which can be irregularly spaced.

Parameters:
  • input_array (numpy array): the array to interpolate

  • x (numpy array): the output coordinates along the x axis

    expressed as (fractional) indices

  • y (numpy array): the output coordinates along the y axis

    expressed as (fractional) indices

  • z (numpy array): the output coordinates along the z axis

    expressed as (fractional) indices

  • order (int): the order of the spline interpolation. Default

    is 0 (linear interpolation). Setting order=1 gives the same behaviour as IDL’s interpolate function with default parameters.

Returns:
Interpolated array with shape (nx, ny, nz), where nx, ny and nz are the lengths of the arrays x, y and z respectively.
c2raytools.smoothing.lanczos_kernel(size, kernel_width)

Generate a 2D Lanczos kernel.

Parameters:
  • size (int): the size of the array
  • kernel_width (int): the width of the kernel
Returns:
The kernel as a (size,size) array
c2raytools.smoothing.smooth_gauss(input_array, sigma)

Smooth the input array with a Gaussian kernel.

Parameters:
  • input_array (numpy array): the array to smooth
  • sigma (float): the width of the kernel
Returns:
The smoothed array. A numpy array with the same dimensions as the input.
c2raytools.smoothing.smooth_lanczos(input_array, kernel_width)

Smooth the input array with a Lanczos kernel.

Parameters:
  • input_array (numpy array): the array to smooth
  • kernel_width (int): the width of the kernel in cells
Returns:
The smoothed array. A numpy array with the same dimensions as the input.
c2raytools.smoothing.smooth_tophat(input_array, tophat_width)

Smooth the input array with a square tophat kernel.

Parameters:
  • input_array (numpy array): the array to smooth
  • tophat_width (int): the width of the kernel in cells
Returns:
The smoothed array. A numpy array with the same dimensions as the input.
c2raytools.smoothing.smooth_with_kernel(input_array, kernel)

Smooth the input array with an arbitrary kernel.

Parameters:
  • input_array (numpy array): the array to smooth

  • kernel (numpy array): the smoothing kernel. Must

    be the same size as the input array

Returns:
The smoothed array. A numpy array with the same dimensions as the input.
c2raytools.smoothing.tophat_kernel(size, tophat_width)

Generate a square tophat kernel

Parameters:
  • size (int): the size of the array
  • tophat_width (int): the size of the tophat kernel
Returns:
The kernel as a (size,size) array
c2raytools.smoothing.tophat_kernel_3d(size)

Generate a 3-dimensional tophat kernel with the specified size

Parameters:
  • size (integer or list-like): the size of

    the tophat kernel along each dimension. If size is an integer, the kernel will be cubic.

Returns:
The normalized kernel

Created on Apr 23, 2015

@author: Hannes Jensen

c2raytools.gaussian_random_field.make_gaussian_random_field(dims, box_dims, power_spectrum, random_seed=None)

Generate a Gaussian random field with the specified power spectrum.

Parameters:
  • dims (tuple): the dimensions of the field in number

    of cells. Can be 2D or 3D.

  • box_dims (float or tuple): the dimensions of the field

    in cMpc.

  • power_spectrum (callable, one parameter): the desired

    spherically-averaged power spectrum of the output. Given as a function of k

  • random_seed (int): the seed for the random number generation

Returns:
The Gaussian random field as a numpy array
c2raytools.gaussian_random_field.make_gaussian_random_field_like_field(input_field, box_dims, random_seed=None)

Generate a Gaussian random field with the same power spectrum as the input field.

Parameters:
  • input_field (numpy array): The field to take the power spectrum

    from

  • box_dims (float or tuple): The dimensions of the input_field in cMpc

  • random_seed (int): the seed for the random number generation

Returns:
A Gaussian random field with the same dimensions and power spectrum as the input field

Note

This function is experimental. There are often interpolation issues at low k. Use with care.