UVPSpec Class

UVPSpec objects contain power spectra and associated metadata.

class hera_pspec.UVPSpec[source]

An object for storing power spectra generated by hera_pspec and a file-format for its data and meta-data.

antnums_to_bl(antnums)[source]

Convert tuple of antenna numbers to baseline integer.

Parameters:antnums (tuple) – Tuple containing integer antenna numbers for a baseline. Ex. (ant1, ant2)
Returns:bl – Baseline integer.
Return type:i6 integer
antnums_to_blpair(antnums)[source]

Convert nested tuple of antenna numbers to baseline-pair integer.

Parameters:antnums (tuple) – nested tuple containing integer antenna numbers for a baseline-pair. Ex. ((ant1, ant2), (ant3, ant4))
Returns:blpair – baseline-pair integer
Return type:i12 integer
average_spectra(blpair_groups=None, time_avg=False, blpair_weights=None, inplace=True)[source]

Average power spectra across the baseline-pair-time axis, weighted by each spectrum’s integration time.

This is an “incoherent” average, in the sense that this averages power spectra, rather than visibility data. The ‘nsample_array’ and ‘integration_array’ will be updated to reflect the averaging.

In the case of averaging across baseline pairs, the resultant averaged spectrum is assigned to the zeroth blpair in the group. In the case of time averaging, the time and LST arrays are assigned to the mean of the averaging window.

Note that this is designed to be separate from spherical binning in k: here we are not connecting k_perp modes to k_para modes. However, if blpairs holds groups of iso baseline separation, then this is equivalent to cylindrical binning in 3D k-space.

If you want help constructing baseline-pair groups from baseline groups, see self.get_blpair_groups_from_bl_groups.

Parameters:

blpair_groups (list of baseline-pair groups) – List of list of tuples or integers. All power spectra in a baseline-pair group are averaged together. If a baseline-pair exists in more than one group, a warning is raised. Examples:

blpair_groups = [ [((1, 2), (1, 2)), ((2, 3), (2, 3))],
                  [((4, 6), (4, 6))]]

blpair_groups = [ [1002001002, 2003002003], [4006004006] ]
time_avg : bool, optional
If True, average power spectra across the time axis. Default: False.
blpair_weights : list of weights (float or int), optional
Relative weight of each baseline-pair when performing the average. This is useful for bootstrapping. This should have the same shape as blpair_groups if specified. The weights are automatically normalized within each baseline-pair group. Default: None (all baseline pairs have unity weights).
inplace : bool, optional
If True, edit data in self, else make a copy and return. Default: True.

Notes

Currently, every baseline-pair in a blpair group must have the same Ntimes, unless time_avg=True. Future versions may support baseline-pair averaging of heterogeneous time arrays. This includes the scenario of repeated blpairs (e.g. in bootstrapping), which will return multiple copies of their time_array.

bl_to_antnums(bl)[source]

Convert baseline (antenna-pair) integer to nested tuple of antenna numbers.

Parameters:bl (i6 int) – baseline integer ID
Returns:antnums – tuple containing baseline antenna numbers. Ex. (ant1, ant2)
Return type:tuple
blpair_to_antnums(blpair)[source]

Convert baseline-pair integer to nested tuple of antenna numbers.

Parameters:blpair (i12 int) – baseline-pair integer ID
Returns:antnums – Nested tuple containing baseline-pair antenna numbers. Ex. ((ant1, ant2), (ant3, ant4))
Return type:tuple
blpair_to_indices(blpair)[source]

Convert a baseline-pair nested tuple ((ant1, ant2), (ant3, ant4)) or a baseline-pair integer into indices to index the blpairts axis of data_array.

Parameters:blpair (tuples or ints) – Nested tuple or blpair i12 integer, Ex. ((1, 2), (3, 4)), or list of blpairs.
check(just_meta=False)[source]

Run checks to make sure metadata and data arrays are properly defined and in the right format. Raises AssertionErrors if checks fail.

Parameters:just_meta (bool, optional) – If True, only check metadata. Default: False.
compute_scalar(spw, pol, num_steps=1000, little_h=True, noise_scalar=False)[source]

Compute power spectrum normalization scalar given an adopted cosmology and a beam model. See pspecbeam.PSpecBeamBase.compute_pspec_scalar for details.

Parameters:
  • spw (integer) – Spectral window selection.
  • pol (str or int) – Which polarization to calculate the scalar for.
  • num_steps (int, optional) – Number of integration bins along frequency in computing scalar. Default: 1000.
  • little_h (bool, optional) – Whether to use h^-1 units or not. Default: True.
  • noise_scalar (bool, optional) – If True calculate noise pspec scalar, else calculate normal pspec scalar. See pspecbeam.py for difference between normal scalar and noise scalar. Default: False.
Returns:

scalar – Power spectrum normalization scalar.

Return type:

float

convert_to_deltasq(little_h=True, inplace=True)[source]

Convert from P(k) to Delta^2(k) by multiplying by k^3 / (2pi^2).

The units of the output is therefore the current units (self.units) times h^3 Mpc^-3, where the h^3 is only included if little_h == True.

Parameters:inplace (boolean) – If True edit and overwrite arrays in self, else make a copy of self and return.
fold_spectra()[source]

Average bandpowers from matching positive and negative delay bins onto a purely positive delay axis. Negative delay bins are still populated, but are filled with zero power. This is an in-place operation.

Will only work if self.folded == False, i.e. data is currently unfolded across negative and positive delay. Because this averages the data, the nsample array is multiplied by a factor of 2.

WARNING: This operation cannot be undone.

generate_noise_spectra(spw, pol, Tsys, blpairs=None, little_h=True, form='Pk', num_steps=2000, real=True)[source]

Generate the expected 1-sigma noise power spectrum given a selection of spectral window, system temp., and polarization. This estimate is constructed as:

P_N = scalar * (Tsys * 1e3)^2 / (integration_time) / sqrt(Nincoherent)
[mK^2 h^-3 Mpc^3]

where scalar is the cosmological and beam scalar (i.e. X2Y * Omega_eff) calculated from pspecbeam with noise_scalar = True, integration_time is in seconds and comes from self.integration_array and Nincoherent is the number of incoherent averaging samples and comes from self.nsample_array.

If the polarization specified is a pseudo Stokes pol (I, Q, U or V) then an extra factor of 2 is divided. If form == ‘DelSq’ then a factor of k^3 / (2pi^2) is multiplied. If real is True, a factor of sqrt(2) is divided to account for discarding imaginary noise component.

For more details, see hera_pspec.noise.Sensitivity.calc_P_N, and Cheng et al. 2018.

The default units of P_N are mK^2 h^-3 Mpc^3 (if little_h=True and form=’Pk’), but is different if those kwargs are altered.

Parameters:
  • spw (int) – Spectral window index to generate noise curve for.
  • pol (str or int) – Polarization selection in form of str (e.g. ‘I’ or ‘Q’ or ‘xx’) or int (e.g. -5 or -6).
  • Tsys (float) – System temperature in Kelvin.
  • blpairs (list) – List of unique blair tuples or i12 integers to calculate noise spectrum for default is to calculate for baseline pairs.
  • little_h (bool, optional) – Whether to have cosmological length units be h^-1 Mpc or Mpc Default: h^-1 Mpc
  • form (str, optional) – Form of power spectra, whetehr P(k) or Delta^2(k). Valid options are ‘Pk’ or ‘DelSq’. Default: ‘Pk’.
  • num_steps (int, optional) – Number of frequency bins to use in integrating power spectrum scalar in pspecbeam. Default: 2000.
  • real (bool, optional) – If True assumes the real component of complex power spectrum is used, and will divide P_N by an extra sqrt(2). Otherwise, assume power spectra are complex and keep P_N as is. Default: True.
Returns:

P_N – Dictionary containing blpair integers as keys and float ndarrays of noise power spectra as values, with ndarrays having shape (Ntimes, Ndlys).

Return type:

dict

get_ENU_bl_vecs()[source]

Return baseline vector array in TOPO (ENU) frame in meters, with matched ordering of self.bl_vecs.

Returns:blvecs – Array of baseline vectors.
Return type:array_like
get_blpair_groups_from_bl_groups(blgroups, only_pairs_in_bls=False)[source]

Get baseline pair matches from a list of baseline groups.

Parameters:
  • blgroups (list) – List of baseline groups, which themselves are lists of baseline tuples or baseline i6 integers. Ex: [ [(1, 2), (2, 3), (3, 4)], [(1, 4), (2, 5)] ]
  • only_pairs_in_bls (bool, optional) – If True, select only baseline-pairs whose first _and_ second baseline are both found in each baseline group. Default: False.
Returns:

blpair_groups – List of blpair groups, which themselves are lists of blpair integers.

Return type:

list

get_blpair_seps()[source]

For each baseline-pair, get the average baseline separation in ENU frame in meters.

Returns:blp_avg_sep – shape=(Nblpairts,)
Return type:float ndarray
get_data(key, *args)[source]

Slice into data_array with a specified data key in the format

(spw, ((ant1, ant2), (ant3, ant4)), pol)

or

(spw, blpair-integer, pol)

where spw is the spectral window integer, ant1 etc. are integers, and pol is either a polarization string (ex. ‘XX’) or integer (ex. -5).

Parameters:key (tuple) – Baseline-pair key
Returns:data – Shape (Ntimes, Ndlys)
Return type:complex ndarray
get_dlys(key)[source]

Get array of delays given a spectral window selection.

Parameters:key (int, or tuple with integer) – Spectral window selection
Returns:dlys
Return type:float ndarray, contains delays in nanosec of pspectra given spw
get_integrations(key, *args)[source]

Slice into integration_array with a specified data key in the format:

(spw, ((ant1, ant2), (ant3, ant4)), pol)

or:

(spw, blpair-integer, pol)

where spw is the spectral window integer, ant1 etc. are integers, and pol is either a polarization string (ex. ‘XX’) or integer (ex. -5).

Parameters:key (tuple) – Baseline-pair key
Returns:data – Has shape (Ntimes,)
Return type:float ndarray
get_kparas(spw, little_h=True)[source]

Get radial (parallel) cosmological wavevectors for power spectra given an adopted cosmology and a spw selection.

Parameters:
  • spw (int, choice of spectral window)
  • little_h (boolean, optional) – Whether to have cosmological length units be h^-1 Mpc or Mpc Default: h^-1 Mpc
  • Returns (k_perp, k_para)
  • ——-
  • k_para (float ndarray, radial wave-vectors, shape=(Ndlys given spw))
get_kperps(spw, little_h=True)[source]

Get transverse (perpendicular) cosmological wavevector for each baseline-pair given an adopted cosmology and a spw selection.

Parameters:
  • spw (int, choice of spectral window)
  • little_h (boolean, optional) – Whether to have cosmological length units be h^-1 Mpc or Mpc Default: h^-1 Mpc
Returns:

k_perp

Return type:

float ndarray, transverse wave-vectors, shape=(Nblpairs,)

get_nsamples(key, *args)[source]

Slice into nsample_array with a specified data key in the format

(spw, ((ant1, ant2), (ant3, ant4)), pol)

or

(spw, blpair-integer, pol)

where spw is the spectral window integer, ant1 etc. are integers, and pol is either a polarization string (ex. ‘XX’) or integer (ex. -5).

Parameters:key (tuple, baseline-pair key)
Returns:data
Return type:float ndarray with shape (Ntimes,)
get_wgts(key, *args)[source]

Slice into wgt_array with a specified data key in the format

(spw, ((ant1, ant2), (ant3, ant4)), pol)

or

(spw, blpair-integer, pol)

where spw is the spectral window integer, ant1 etc. are integers, and pol is either a polarization string (ex. ‘XX’) or integer (ex. -5).

Parameters:key (tuple, baseline-pair key)
Returns:wgts – Has shape (2, Ntimes, Ndlys), where the zeroth axis holds [wgt_1, wgt_2] in that order.
Return type:float ndarray
key_to_indices(key, *args)[source]

Convert a data key into relevant slice arrays. A data key takes the form

(spw_integer, ((ant1, ant2), (ant3, ant4)), pol_string)

or

(spw, blpair-integer, pol)

where spw is the spectral window integer, ant1 etc. are integers, and pol is either a polarization string (ex. ‘XX’) or integer (ex. -5).

One can also expand this key into the kwarg slots, such that key=spw, key2=blpair, and key3=pol.

The key can also be a dictionary in the form:

key = {
  'spw' : spw_integer,
  'blpair' : ((ant1, ant2), (ant3, ant4))
  'pol' : pol_string
  }

and it will parse the dictionary for you.

Parameters:key (tuple) – Baseline-pair key.
Returns:
  • spw (int) – Spectral window index.
  • blpairts (list) – List of integers to apply along blpairts axis.
  • pol (int) – Polarization index.
pol_to_indices(pol)[source]

Map a polarization integer or str to its index in pol_array

Parameters:pol (str or int) – Polarization string (ex. ‘XX’) or integer (ex. -5), or a list of strs or ints.
Returns:indices – Index of pol in pol_array.
Return type:int
read_from_group(grp, just_meta=False, spws=None, bls=None, blpairs=None, times=None, pols=None, only_pairs_in_bls=False)[source]

Clear current UVPSpec object and load in data from specified HDF5 group.

Parameters:
  • grp (HDF5 group) – HDF5 group to load data from.
  • just_meta (bool, optional) – If True, read-in metadata but ignore data, wgts and integration arrays. Default: False.
  • spws (list of tuple, optional) – List of spectral window integers to select. Default: None (loads all channels).
  • bls (list of i6 baseline integers or baseline tuples) – Select all baseline-pairs whose first _or_ second baseline are in the list. This changes if only_pairs_in_bls == True. Example tuple: (2, 3). Default: None (loads all bl pairs).
  • blpairs (list of baseline-pair tuples or integers) – List of baseline pairs to keep. If bls is also fed, this list is concatenated onto the baseline-pair list constructed from from the bls selection.
  • times (float ndarray) – Times from the time_avg_array to keep.
  • pols (list of str or int) – List of polarization strings or integers to keep. See pyuvdata.utils.polstr2num for acceptable options.
  • only_pairs_in_bls (bool, optional) – If True, keep only baseline-pairs whose first _and_ second baseline are found in the ‘bls’ list. Default: False.
read_hdf5(filepath, just_meta=False, spws=None, bls=None, blpairs=None, times=None, pols=None, only_pairs_in_bls=False)[source]

Clear current UVPSpec object and load in data from an HDF5 file.

Parameters:
  • filepath (str) – Path to HDF5 file to read from.
  • just_meta (bool, optional) – If True, read-in metadata but ignore data, wgts and integration arrays. Default: False.
  • spws (list of tuple, optional) – List of spectral window integers to select. Default: None (loads all channels).
  • bls (list of i6 baseline integers or baseline tuples) – Select all baseline-pairs whose first _or_ second baseline are in the list. This changes if only_pairs_in_bls == True. Example tuple: (2, 3). Default: None (loads all bl pairs).
  • blpairs (list of baseline-pair tuples or integers) – List of baseline pairs to keep. If bls is also fed, this list is concatenated onto the baseline-pair list constructed from from the bls selection.
  • times (float ndarray) – Times from the time_avg_array to keep.
  • pols (list of str or int) – List of polarization strings or integers to keep. See pyuvdata.utils.polstr2num for acceptable options.
  • only_pairs_in_bls (bool, optional) – If True, keep only baseline-pairs whose first _and_ second baseline are found in the ‘bls’ list. Default: False.
select(spws=None, bls=None, only_pairs_in_bls=True, blpairs=None, times=None, pols=None, inplace=True)[source]

Select function for selecting out certain slices of the data.

Parameters:
  • spws (list, optional) – List of spectral window integers to select. If None, all will be selected. Default: None.
  • bls (list of i6 baseline integers or baseline tuples, optional) – Example: (2,3). Select all baseline-pairs whose first _or_ second baseline are in bls list. This changes if only_pairs_in_bls == True. If None, all baselines are selected (subject to other options). Default: None.
  • only_pairs_in_bls (bool, optional) – If True, keep only baseline-pairs whose first _and_ second baseline are found in bls list. Default: True.
  • blpairs (list of baseline-pair tuples or integers, optional) – List of baseline-pairs to keep. If bls is also fed, this list is concatenated onto the baseline-pair list constructed from the bls selection. If None, all valid baseline pairs are selected. Default: None.
  • times (array_like, optional) – Float ndarray of times from the time_avg_array to select. If None, all times are kept. Default: None.
  • pols (list of str or int, optional) – List of polarizations to keep. See pyuvdata.utils.polstr2num for acceptable options. If None, all polarizations are kept. Default: None.
  • inplace (bool, optional) – If True, edit and overwrite arrays in self, else make a copy of self and return. Default: True.
Returns:

uvp – If inplace=False, return a new UVPSpec object containing only the selected data.

Return type:

UVPSpec, optional

set_cosmology(new_cosmo, overwrite=False, new_beam=None, verbose=True)[source]

Set the cosmology for this UVPSpec object by passing an instance of conversions.Cosmo_Conversions and assigning it to self.cosmo, in addition to re-computing power spectrum normalization quantities like self.scalar_array. Because this will attempt to recompute the scalar_array, the beam-related metadata (OmegaP, OmegaPP and beam_freqs) must exist in self.

If they do not, or if you’d like to overwrite them with a new beam, you can pass a UVBeam object or path to a beam file via the new_beam kwarg.

If self.cosmo already exists then you are attempting to overwrite the currently adopted cosmology, which will only proceed if overwrite is True.

If overwrite == True, then this module will overwrite self.cosmo. It will also recompute the power spectrum normalization scalar (using pspecbeam) and will overwrite the values in self.scalar_array. It will then propagate these re-normalization changes into the data_array by multiplying the data by (new_scalar_array / old_scalar_array).

Parameters:
  • new_cosmo (Cosmo_Conversions or dict) – Cosmo_Conversions instance or cosmological parameter dictionary. The new cosmology you want to adopt for this UVPSpec object.
  • overwrite (bool, optional) – If True, overwrite self.cosmo if it already exists. Default: False.
  • new_beam (PSpecBeamUV or str) – pspecbeam.PSpecBeamUV object or path to beam file. The new beam you want to adopt for this UVPSpec object.
  • verbose (bool, optional) – If True, report feedback to stdout. Default: True.
spw_to_indices(spw)[source]

Convert a spectral window integer into a list of indices to index into the spwdlys axis of dly_array and/or freq_array.

If self.folded == False, return indices for both positive and negative delay bins, else if self.folded == True, returns indices for only positive delay bins.

Parameters:spw (int) – Spectral window index or list of indices.
time_to_indices(time, blpairs=None)[source]

Convert a time [Julian Date] from self.time_avg_array to an array that indexes the elements at which it appears in that array. Can optionally take a blpair selection to further select the indices at which both that time and blpair are satisfied.

Parameters:
  • time (float) – Julian Date time from self.time_avg_array, Ex: 2458042.12242
  • blpairs (tuple or int, optional) – List of blpair tuples or integers that further selects the elements at which both time and blpairs are satisfied.
Returns:

indices – Contains indices at which selection is satisfied along the blpairts axis.

Return type:

integer ndarray

units

Return power spectrum units. See self.vis_units and self.norm_units.

write_hdf5(filepath, overwrite=False, run_check=True)[source]

Write a UVPSpec object to HDF5 file.

Parameters:
  • filepath (str) – Filepath for output file.
  • overwrite (bool, optional) – Whether to overwrite output file if it exists. Default: False.
  • run_check (bool, optional) – Run UVPSpec validity check before writing to file. Default: True.
write_to_group(group, run_check=True)[source]

Write UVPSpec data into an HDF5 group.

Parameters:
  • group (HDF5 group) – The handle of the HDF5 group that the UVPSpec object should be written into.
  • run_check (bool, optional) – Whether to run a validity check on the UVPSpec object before writing it to the HDF5 group. Default: True.