PSpecData Class

PSpecData takes UVData objects and calculates delay power spectra from them.

class hera_pspec.PSpecData(dsets=[], wgts=None, dsets_std=None, labels=None, beam=None)[source]
C_model(key, model='empirical')[source]

Return a covariance model having specified a key and model type.

Parameters:
  • key (tuple) – Tuple containing indices of dataset and baselines. The first item specifies the index (ID) of a dataset in the collection, while subsequent indices specify the baseline index, in _key2inds format.
  • model (string, optional) – Type of covariance model to calculate, if not cached. options=[‘empirical’]
Returns:

C – Covariance model for the specified key.

Return type:

array-like

I(key)[source]

Return identity covariance matrix.

Parameters:key (tuple) – Tuple containing indices of dataset and baselines. The first item specifies the index (ID) of a dataset in the collection, while subsequent indices specify the baseline index, in _key2inds format.
Returns:I – Identity covariance matrix, dimension (Nfreqs, Nfreqs).
Return type:array_like
Jy_to_mK(beam=None)[source]

Convert internal datasets from a Jy-scale to mK scale using a primary beam model if available. Note that if you intend to rephase_to_dset(), Jy to mK conversion must be done after that step.

Parameters:beam (PSpecBeam object) – Beam object.
R(key)[source]

Return the data-weighting matrix R, which is a product of data covariance matrix (I or C^-1), diagonal flag matrix (Y) and diagonal tapering matrix (T):

R = sqrt(T^t) sqrt(Y^t) K sqrt(Y) sqrt(T)

where T is a diagonal matrix holding the taper and Y is a diagonal matrix holding flag weights. The K matrix comes from either I or iC depending on self.data_weighting, T is informed by self.taper and Y is taken from self.Y().

Parameters:key (tuple) – Tuple containing indices of dataset and baselines. The first item specifies the index (ID) of a dataset in the collection, while subsequent indices specify the baseline index, in _key2inds format.
Y(key)[source]

Return the weighting (diagonal) matrix, Y. This matrix is calculated by taking the logical AND of flags across all times given the dset-baseline-pol specification in ‘key’, converted into a float, and inserted along the diagonal of an spw_Nfreqs x spw_Nfreqs matrix.

The logical AND step implies that all time-dependent flagging patterns are automatically broadcasted across all times. This broadcasting follows the principle that, for each freq channel, if at least a single time is unflagged, then the channel is treated as unflagged for all times. Power spectra from certain times, however, can be given zero weight by setting the nsample array to be zero at those times (see self.broadcast_dset_flags).

Parameters:key (tuple) – Tuple containing indices of dataset and baselines. The first item specifies the index (ID) of a dataset in the collection, while subsequent indices specify the baseline index, in _key2inds format.
Returns:Y – spw_Nfreqs x spw_Nfreqs diagonal matrix holding AND of flags across all times for each freq channel.
Return type:array_like
add(dsets, wgts, labels=None, dsets_std=None)[source]

Add a dataset to the collection in this PSpecData object.

Parameters:
  • dsets (UVData or list or dict) – UVData object or list of UVData objects containing data to add to the collection.
  • wgts (UVData or list or dict) – UVData object or list of UVData objects containing weights to add to the collection. Must be the same length as dsets. If a weight is set to None, the flags of the corresponding dset are used.
  • labels (list of str) – An ordered list of names/labels for each dataset, if dsets was specified as a list. If dsets was specified as a dict, the keys of that dict will be used instead.
  • dsets_std (UVData or list or dict) – Optional UVData object or list of UVData objects containing the standard deviations (real and imaginary) of data to add to the collection. If dsets is a dict, will assume dsets_std is a dict and if dsets is a list, will assume dsets_std is a list.
broadcast_dset_flags(spw_ranges=None, time_thresh=0.2, unflag=False)[source]

For each dataset in self.dset, update the flag_array such that the flagging patterns are time-independent for each baseline given a selection for spectral windows.

For each frequency pixel in a selected spw, if the fraction of flagged times exceeds time_thresh, then all times are flagged. If it does not, the specific integrations which hold flags in the spw are flagged across all frequencies in the spw.

Additionally, one can also unflag the flag_array entirely if desired.

Note: although technically allowed, this function may give unexpected results if multiple spectral windows in spw_ranges have frequency overlap.

Note: it is generally not recommended to set time_thresh > 0.5, which could lead to substantial amounts of data being flagged.

Parameters:
  • spw_ranges (list of tuples) – list of len-2 spectral window tuples, specifying the start (inclusive) and stop (exclusive) index of the frequency array for each spw. Default is to use the whole band
  • time_thresh (float) – Fractional threshold of flagged pixels across time needed to flag all times per freq channel. It is not recommend to set this greater than 0.5
  • unflag (bool) – If True, unflag all data in the spectral window.
check_key_in_dset(key, dset_ind)[source]

Check ‘key’ exists in the UVData object self.dsets[dset_ind]

Parameters:
  • key (tuple) – if length 1: assumed to be polarization number or string elif length 2: assumed to be antenna-number tuple (ant1, ant2) elif length 3: assumed ot be antenna-number-polarization tuple (ant1, ant2, pol)
  • dset_ind (int, the index of the dataset to-be-checked)
Returns:

exists – True if the key exists, False otherwise

Return type:

bool

clear_cache(keys=None)[source]

Clear stored matrix data (or some subset of it).

Parameters:keys (list of tuples, optional) – List of keys to remove from matrix cache. If None, all keys will be removed. Default: None.
cov_p_hat(M, q_cov)[source]

Covariance estimate between two different band powers p_alpha and p_beta given by M_{alpha i} M^*_{beta,j} C_q^{ij} where C_q^{ij} is the q-covariance

Parameters:
  • M (array_like) – Normalization matrix, M.
  • q_cov (array_like) – covariance between bandpowers in q_alpha and q_beta
cov_q_hat(key1, key2, time_indices=None)[source]

Compute the un-normalized covariance matrix for q_hat for a given pair of visibility vectors. Returns the following matrix:

Cov(hat{q}_a,hat{q}_b)

!!!Only supports covariance between same power-spectrum estimates!!! (covariance between pair of baselines with the same pair of baselines) !!!Assumes that both baselines used in power-spectrum estimate !!!have independent noise relizations!!!

Parameters:
  • key1, key2 (tuples or lists of tuples) – Tuples containing the indices of the dataset and baselines for the two input datavectors. If a list of tuples is provided, the baselines in the list will be combined with inverse noise weights.
  • time_indices (list of indices of times to include or just a single time.)
  • default is None -> compute covariance for all times.
Returns:

cov_q_hat

Return type:

matrix with covariances between un-normalized band powers

cross_covar_model(key1, key2, model='empirical', conj_1=False, conj_2=True)[source]

Return a covariance model having specified a key and model type.

Parameters:
  • key1, key2 (tuples) – Tuples containing indices of dataset and baselines. The first item specifies the index (ID) of a dataset in the collection, while subsequent indices specify the baseline index, in _key2inds format.
  • model (string, optional) – Type of covariance model to calculate, if not cached. options=[‘empirical’]
  • conj_1 (boolean, optional) – Whether to conjugate first copy of data in covar or not. Default: False
  • conj_2 (boolean, optional) – Whether to conjugate second copy of data in covar or not. Default: True
Returns:

cross_covar – Cross covariance model for the specified key.

Return type:

array-like, spw_Nfreqs x spw_Nfreqs

delays()[source]

Return an array of delays, tau, corresponding to the bins of the delay power spectrum output by pspec() using self.spw_range to specify the spectral window.

Returns:delays – Delays, tau. Units: ns.
Return type:array_like
dset_idx(dset)[source]

Return the index of a dataset, regardless of whether it was specified as an integer of a string.

Parameters:dset (int or str) – Index or name of a dataset belonging to this PSpecData object.
Returns:dset_idx – Index of dataset.
Return type:int
dx(key)[source]

Get standard deviation of data for given dataset and baseline as pecified in standard key format.

Parameters:key (tuple) – Tuple containing datset ID and baseline index. The first element of the tuple is the dataset index (or label), and the subsequent elements are the baseline ID.
Returns:dx – Array of std data from the requested UVData dataset and baseline.
Return type:array_like
get_G(key1, key2)[source]

Calculates

G_ab = (1/2) Tr[R_1 Q^alt_a R_2 Q^alt_b],

which is needed for normalizing the power spectrum (see HERA memo #44).

Note that in the limit that R_1 = R_2 = C^-1, this reduces to the Fisher matrix

F_ab = 1/2 Tr [C^-1 Q^alt_a C^-1 Q^alt_b] (arXiv:1502.06016, Eq. 17)
Parameters:key1, key2 (tuples or lists of tuples) – Tuples containing indices of dataset and baselines for the two input datavectors. If a list of tuples is provided, the baselines in the list will be combined with inverse noise weights.
Returns:G – Fisher matrix, with dimensions (Nfreqs, Nfreqs).
Return type:array_like, complex
get_H(key1, key2, sampling=False)[source]

Calculates the response matrix H of the unnormalized band powers q to the true band powers p, i.e.,

<q_a> = sum_b H_{ab} p_b

This is given by

H_ab = (1/2) Tr[R_1 Q_a^alt R_2 Q_b]

(See HERA memo #44). As currently implemented, this approximates the primary beam as frequency independent. Under this approximation, the our H_ab is defined using the equation above except we have Q^tapered rather than Q_b, where

overline{Q}^{tapered,beta} = e^{i 2pi eta_beta (nu_i - nu_j)} gamma(nu_i) gamma(nu_j)

where gamma is the tapering function. Again, see HERA memo #44 for details.

The sampling option determines whether one is assuming that the output points are integrals over k bins or samples at specific k values. The effect is to add a dampening of widely separated frequency correlations with the addition of the term

sinc(pi Delta eta (nu_i - nu_j))

Note that numpy uses the engineering definition of sinc, where sinc(x) = sin(pi x) / (pi x), whereas in the line above and in all of our documentation, we use the physicist definition where sinc(x) = sin(x) / x

Note that in the limit that R_1 = R_2 = C^-1 and Q_a is used instead of Q_a^alt, this reduces to the Fisher matrix

F_ab = 1/2 Tr [C^-1 Q_a C^-1 Q_b] (arXiv:1502.06016, Eq. 17)

This function uses the state of self.taper in constructing H. See PSpecData.pspec for details.

Parameters:
  • key1, key2 (tuples or lists of tuples) – Tuples containing indices of dataset and baselines for the two input datavectors. If a list of tuples is provided, the baselines in the list will be combined with inverse noise weights.
  • sampling (boolean, optional) – Whether to sample the power spectrum or to assume integrated bands over wide delay bins. Default: False
Returns:

H – Dimensions (Nfreqs, Nfreqs).

Return type:

array_like, complex

get_MW(G, H, mode='I', band_covar=None)[source]

Construct the normalization matrix M and window function matrix W for the power spectrum estimator. These are defined through Eqs. 14-16 of arXiv:1502.06016:

hat{p} = M hat{q} <hat{p}> = W p W = M H,

where p is the true band power and H is the response matrix (defined above in get_H) of unnormalized bandpowers to normed bandpowers.

Several choices for M are supported:
‘I’: Set M to be diagonal (e.g. HERA Memo #44) ‘H^-1’: Set M = H^-1, the (pseudo)inverse response matrix. ‘V^-1/2’: Set M = V^-1/2, the root-inverse response matrix (using SVD).
These choices will be supported very soon:
‘L^-1’: Set M = L^-1, Cholesky decomposition.

As written, the window functions will not be correclty normalized; it needs to be adjusted by the pspec scalar for it to be approximately correctly normalized. If the beam is being provided, this will be done in the pspec function.

Parameters:
  • G (array_like) – Denominator matrix for the bandpowers, with dimensions (Nfreqs, Nfreqs).
  • H (array_like) – Response matrix for the bandpowers, with dimensions (Nfreqs, Nfreqs).
  • mode (str, optional) – Definition to use for M. Must be one of the options listed above. Default: ‘I’.
  • band_covar (array_like, optional) – Covariance matrix of the unnormalized bandpowers (i.e., q). Used only if requesting the V^-1/2 normalization. Use get_unnormed_V to get the covariance to put in here, or provide your own array. Default: None
Returns:

  • M (array_like) – Normalization matrix, M. (If G was passed in as a dict, a dict of array_like will be returned.)
  • W (array_like) – Window function matrix, W. (If G was passed in as a dict, a dict of array_like will be returned.)

get_Q(mode, pol=False)[source]

Computes Q_alpha(i,j), which is the response of the data covariance to the bandpower (dC/dP_alpha). This includes contributions from primary beam.

Parameters:
  • mode (int) – Central wavenumber (index) of the bandpower, p_alpha.
  • pol (str/int/bool, optional) – Which beam polarization to use. If the specified polarization doesn’t exist, a uniform isotropic beam (with integral 4pi for all frequencies) is assumed. Default: False (uniform beam).
Returns:

Q – Response matrix for bandpower p_alpha.

Return type:

array_like

get_Q_alt(mode, allow_fft=True)[source]

Response of the covariance to a given bandpower, dC / dp_alpha, EXCEPT without the primary beam factors. This is Q_alt as defined in HERA memo #44, so it’s not dC / dp_alpha, strictly, but is just the part that does the Fourier transforms.

Assumes that Q will operate on a visibility vector in frequency space. In the limit that self.spw_Ndlys equals self.spw_Nfreqs, this will produce a matrix Q that performs a two-sided FFT and extracts a particular Fourier mode.

(Computing x^t Q y is equivalent to Fourier transforming x and y separately, extracting one element of the Fourier transformed vectors, and then multiplying them.)

When self.spw_Ndlys < self.spw_Nfreqs, the effect is similar except the delay bins need not be in the locations usually mandated by the FFT algorithm.

Parameters:
  • mode (int) – Central wavenumber (index) of the bandpower, p_alpha.
  • allow_fft (boolean, optional) – If set to True, allows a shortcut FFT method when the number of delay bins equals the number of delay channels. Default: True
Returns:

Q – Response matrix for bandpower p_alpha.

Return type:

array_like

get_unnormed_E(key1, key2)[source]

Calculates a series of unnormalized E matrices, such that

q_a = x_1^* E^{12,a} x_2

so that

E^{12,a} = (1/2) R_1 Q^a R_2.

In principle, this could be used to actually estimate q-hat. In other words, we could call this function to get E, and then sandwich it with two data vectors to get q_a. However, this should be slower than the implementation in q_hat. So this should only be used as a helper method for methods such as get_unnormed_V.

Note for the future: There may be advantages to doing the matrix multiplications separately

Parameters:key1, key2 (tuples or lists of tuples) – Tuples containing indices of dataset and baselines for the two input datavectors. If a list of tuples is provided, the baselines in the list will be combined with inverse noise weights.
Returns:E – Set of E matrices, with dimensions (Ndlys, Nfreqs, Nfreqs).
Return type:array_like, complex
get_unnormed_V(key1, key2, model='empirical')[source]

Calculates the covariance matrix for unnormed bandpowers (i.e., the q vectors). If the data were real and x_1 = x_2, the expression would be

\[V_ab = 2 tr(C E_a C E_b), where E_a = (1/2) R Q^a R\]

When the data are complex, the expression becomes considerably more complicated. Define

\[E^{12,a} = (1/2) R_1 Q^a R_2 C^1 = <x1 x1^dagger> - <x1><x1^dagger> C^2 = <x2 x2^dagger> - <x2><x2^dagger> P^{12} = <x1 x2> - <x1><x2> S^{12} = <x1^* x2^*> - <x1^*> <x2^*>\]

Then

\[V_ab = tr(E^{12,a} C^2 E^{21,b} C^1) + tr(E^{12,a} P^{21} E^{12,b *} S^{21})\]

Note that

\[E^{12,a}_{ij}.conj = E^{21,a}_{ji}\]

This function estimates C^1, C^2, P^{12}, and S^{12} empirically by default. (So while the pointy brackets <…> should in principle be ensemble averages, in practice the code performs averages in time.)

Empirical covariance estimates are in principle a little risky, as they can potentially induce signal loss. This is probably ok if we are just looking intending to look at V. It is most dangerous when C_emp^-1 is applied to the data. The application of using this to form do a V^-1/2 decorrelation is probably medium risk. But this has yet to be proven, and results coming from V^-1/2 should be interpreted with caution.

Note for future: Although the V matrix should be Hermitian by construction, in practice there are precision issues and the Hermiticity is violated at ~ 1 part in 10^15. (Which is ~the expected roundoff error). If something messes up, it may be worth investigating this more.

Note for the future: If this ends up too slow, Cholesky tricks can be employed to speed up the computation by a factor of a few.

Parameters:
  • key1, key2 (tuples or lists of tuples) – Tuples containing indices of dataset and baselines for the two input datavectors. If a list of tuples is provided, the baselines in the list will be combined with inverse noise weights.
  • model (str, default: ‘empirical’) – How the covariances of the input data should be estimated.
Returns:

V – Bandpower covariance matrix, with dimensions (Nfreqs, Nfreqs).

Return type:

array_like, complex

iC(key, model='empirical')[source]

Return the inverse covariance matrix, C^-1.

Parameters:
  • key (tuple) – Tuple containing indices of dataset and baselines. The first item specifies the index (ID) of a dataset in the collection, while subsequent indices specify the baseline index, in _key2inds format.
  • model (string) – Type of covariance model to calculate, if not cached. options=[‘empirical’]
Returns:

iC – Inverse covariance matrix for specified dataset and baseline.

Return type:

array_like

p_hat(M, q)[source]

Optimal estimate of bandpower p_alpha, defined as p_hat = M q_hat.

Parameters:
  • M (array_like) – Normalization matrix, M.
  • q (array_like) – Unnormalized bandpowers, hat{q}.
Returns:

p_hat – Optimal estimate of bandpower, hat{p}.

Return type:

array_like

parse_blkey(key)[source]

Parse a dataset + baseline key in the form (dataset, baseline, [pol]) into a dataset index and a baseline(pol) key, where pol is optional.

Parameters:key (tuple)
Returns:
  • dset (int) – dataset index
  • bl (tuple) – baseline (pol) key
pspec(bls1, bls2, dsets, pols, n_dlys=None, input_data_weight='identity', norm='I', taper='none', sampling=False, little_h=True, spw_ranges=None, baseline_tol=1.0, store_cov=False, verbose=True, exact_norm=False, history='')[source]

Estimate the delay power spectrum from a pair of datasets contained in this object, using the optimal quadratic estimator of arXiv:1502.06016.

In this formulation, the power spectrum is proportional to the visibility data via

P = M data_{LH} E data_{RH}

where E contains the data weighting and FT matrices, M is a normalization matrix, and the two separate datasets are denoted as “left-hand” and “right-hand”.

Each power spectrum is generated by taking a baseline (specified by an antenna-pair and polarization key) from bls1 out of dsets[0] and assigning it as data_LH, and a bl from bls2 out of the dsets[1] and assigning it as data_RH.

If the bl chosen from bls1 is (ant1, ant2) and the bl chosen from bls2 is (ant3, ant4), the “baseline-pair” describing their cross multiplication is ((ant1, ant2), (ant3, ant4)).

Parameters:
  • bls1, bls2 (list) – List of baseline groups, each group being a list of ant-pair tuples.

  • dsets (length-2 tuple or list) – Contains indices of self.dsets to use in forming power spectra, where the first index is for the Left-Hand dataset and second index is used for the Right-Hand dataset (see above).

  • pols (tuple or list of tuple) – Contains polarization pairs to use in forming power spectra e.g. (‘XX’,’XX’) or [(‘XX’,’XX’),(‘pI’,’pI’)] or a list of polarization pairs. Individual strings are also supported, and will be expanded into a matching pair of polarizations, e.g. ‘xx’ becomes (‘xx’, ‘xx’).

    If a primary_beam is specified, only equal-polarization pairs can be cross-correlated, as the beam scalar normalization is only implemented in this case. To obtain unnormalized spectra for pairs of different polarizations, set the primary_beam to None.

  • n_dlys (list of integer, optional) – The number of delay bins to use. The order in the list corresponds to the order in spw_ranges. Default: None, which then sets n_dlys = number of frequencies.

  • input_data_weight (str, optional) – String specifying which weighting matrix to apply to the input data. See the options in the set_R() method for details. Default: ‘identity’.

  • norm (str, optional) – String specifying how to choose the normalization matrix, M. See the ‘mode’ argument of get_MW() for options. Default: ‘I’.

  • taper (str, optional) – Tapering (window) function to apply to the data. Takes the same arguments as aipy.dsp.gen_window(). Default: ‘none’.

  • sampling (boolean, optional) – Whether output pspec values are samples at various delay bins or are integrated bandpowers over delay bins. Default: False

  • little_h (boolean, optional) – Whether to have cosmological length units be h^-1 Mpc or Mpc Default: h^-1 Mpc

  • spw_ranges (list of tuples, optional) – A list of spectral window channel ranges to select within the total bandwidth of the datasets, each of which forms an independent power spectrum estimate. Example: [(220, 320), (650, 775)]. Each tuple should contain a start (inclusive) and stop (exclusive) channel used to index the freq_array of each dataset. The default (None) is to use the entire band provided in each dataset.

  • baseline_tol (float, optional) – Distance tolerance for notion of baseline “redundancy” in meters. Default: 1.0.

  • store_cov (boolean, optional) – If True, calculate an analytic covariance between bandpowers given an input visibility noise model, and store the output in the UVPSpec object.

  • verbose (bool, optional) – If True, print progress, warnings and debugging info to stdout.

  • exact_norm (bool, optional) – If True, estimates power spectrum using Q instead of Q_alt (HERA memo #44). The default options is False. Beware that turning this True would take ~ 7 sec for computing power spectrum for 100 channels per time sample per baseline.

  • history (str, optional) – history string to attach to UVPSpec object

Returns:

uvp – Instance of UVPSpec that holds the output power spectrum data.

Return type:

UVPSpec object

Examples

Example 1: No grouping; i.e. each baseline is its own group, no brackets needed for each bl. If:

A = (1, 2); B = (2, 3); C = (3, 4); D = (4, 5); E = (5, 6); F = (6, 7)

and:

bls1 = [ A, B, C ]
bls2 = [ D, E, F ]

then:

blpairs = [ (A, D), (B, E), (C, F) ]

Example 2: Grouping; blpairs come in lists of blgroups, which are considered “grouped” in OQE. If:

bls1 = [ [A, B], [C, D] ]
bls2 = [ [C, D], [E, F] ]

then:

blpairs = [ [(A, C), (B, D)], [(C, E), (D, F)] ]

Example 3: Mixed grouping; i.e. some blpairs are grouped, others are not. If:

bls1 = [ [A, B], C ]
bls2 = [ [D, E], F ]

then:

blpairs = [ [(A, D), (B, E)], (C, F)]
q_hat(key1, key2, allow_fft=False, exact_norm=False, pol=False)[source]

If exact_norm is False: Construct an unnormalized bandpower, q_hat, from a given pair of visibility vectors. Returns the following quantity:

hat{q}_a = (1/2) conj(x_1) R_1 Q^alt_a R_2 x_2

Note that the R matrix need not be set to C^-1. This is something that is set by the user in the set_R method.

This is related to Equation 13 of arXiv:1502.06016. However, notice that there is a Q^alt_a instead of Q_a. The latter is defined as Q_a equiv dC/dp_a. Since this is the derivative of the covariance with respect to the power spectrum, it contains factors of the primary beam etc. Q^alt_a strips away all of this, leaving only the barebones job of taking a Fourier transform. See HERA memo #44 for details.

This function uses the state of self.data_weighting and self.taper in constructing q_hat. See PSpecData.pspec for details.

If exact_norm is True: It uses get_Q function to return normalized power spectrum (Eq. 14 in HERA memo #44)

Parameters:
  • key1, key2 (tuples or lists of tuples) – Tuples containing indices of dataset and baselines for the two input datavectors for each power spectrum estimate. q_a formed from key1, key2 If a list of tuples is provided, the baselines in the list will be combined with inverse noise weights.
  • allow_fft (bool, optional) – Whether to use a fast FFT summation trick to construct q_hat, or a simpler brute-force matrix multiplication. The FFT method assumes a delta-fn bin in delay space. It also only works if the number of delay bins is equal to the number of frequencies. Default: False.
  • exact_norm (bool, optional) – If True, beam and spectral window factors are taken in the computation of Q_matrix (dC/dp = Q, and not Q_alt) (HERA memo #44, Eq. 11). Q matrix, for each delay mode, is weighted by the integral of beam over theta,phi. Therefore the output power spectra is, by construction, normalized. If True, it returns normalized power spectrum, except for X2Y term. If False, Q_alt is used (HERA memo #44, Eq. 16), and the power spectrum is normalized separately.
  • pol (str/int/bool, optional) – Used only if exact_norm is True. This argument is passed to get_Q to extract the requested beam polarization. Default is the first polarization passed to pspec.
Returns:

q_hat – Unnormalized/normalized bandpowers

Return type:

array_like

rephase_to_dset(dset_index=0, inplace=True)[source]

Rephase visibility data in self.dsets to the LST grid of dset[dset_index] using hera_cal.utils.lst_rephase.

Each integration in all other dsets is phased to the center of the corresponding LST bin (by index) in dset[dset_index].

Will only phase if the dataset’s phase type is ‘drift’. This is because the rephasing algorithm assumes the data is drift-phased when applying the phasor term.

Note that PSpecData.Jy_to_mK() must be run after rephase_to_dset(), if one intends to use the former capability at any point.

Parameters:
  • dset_index (int or str) – Index or label of dataset in self.dset to phase other datasets to.
  • inplace (bool, optional) – If True, edits data in dsets in-memory. Else, makes a copy of dsets, edits data in the copy and returns to user.
Returns:

  • if inplace – return new_dsets
  • else – return None

scalar(polpair, little_h=True, num_steps=2000, beam=None, taper_override='no_override', exact_norm=False)[source]

Computes the scalar function to convert a power spectrum estimate in “telescope units” to cosmological units, using self.spw_range to set spectral window.

See arxiv:1304.4991 and HERA memo #27 for details.

This function uses the state of self.taper in constructing scalar. See PSpecData.pspec for details.

Parameters:
  • polpair (tuple, int, or str) – Which pair of polarizations to compute the beam scalar for, e.g. (‘pI’, ‘pI’) or (‘XX’, ‘YY’). If string, will assume that the specified polarization is to be cross-correlated with itself, e.g. ‘XX’ implies (‘XX’, ‘XX’).
  • little_h (boolean, optional) – Whether to have cosmological length units be h^-1 Mpc or Mpc Default: h^-1 Mpc
  • num_steps (int, optional) – Number of steps to use when interpolating primary beams for numerical integral Default: 10000
  • beam (PSpecBeam object) – Option to use a manually-fed PSpecBeam object instead of using self.primary_beam.
  • taper_override (str, optional) – Option to override the taper chosen in self.taper (does not overwrite self.taper; just applies to this function). Default: no_override
  • exact_norm (boolean, optional) – If True, scalar would just be the X2Y term, as the beam and spectral terms are taken into account while constructing Q matrix.
Returns:

scalar – [int dnu (Omega_PP / Omega_P^2) ( B_PP / B_P^2 ) / (X^2 Y)]^-1 in h^-3 Mpc^3 or Mpc^3.

Return type:

float

scalar_delay_adjustment(key1=None, key2=None, sampling=False, Gv=None, Hv=None)[source]

Computes an adjustment factor for the pspec scalar that is needed when the number of delay bins is not equal to the number of frequency channels.

This adjustment is necessary because sum_gamma tr[Q^alt_alpha Q^alt_gamma] = N_freq**2 is something that is true only when N_freqs = N_dlys.

In general, the result is still independent of alpha, but is no longer given by N_freq**2. (Nor is it just N_dlys**2!)

This function uses the state of self.taper in constructing adjustment. See PSpecData.pspec for details.

Parameters:
  • key1, key2 (tuples or lists of tuples, optional) – Tuples containing indices of dataset and baselines for the two input datavectors. If a list of tuples is provided, the baselines in the list will be combined with inverse noise weights. If Gv and Hv are specified, these arguments will be ignored. Default: None.
  • sampling (boolean, optional) – Whether to sample the power spectrum or to assume integrated bands over wide delay bins. Default: False
  • Gv, Hv (array_like, optional) – If specified, use these arrays instead of calling self.get_G() and self.get_H(). Using precomputed Gv and Hv will speed up this function significantly. Default: None.
Returns:

adjustment

Return type:

float

set_C(cov)[source]

Set the cached covariance matrix to a set of user-provided values.

Parameters:cov (dict) – Dictionary containing new covariance values for given datasets and baselines. Keys of the dictionary are tuples, with the first item being the ID (index) of the dataset, and subsequent items being the baseline indices.
set_Ndlys(ndlys=None)[source]

Set the number of delay bins used.

Parameters:ndlys (integer) – Number of delay bins. Default: None, sets number of delay bins equal to the number of frequency channels
set_R(d)[source]

Set the data-weighting matrix for a given dataset and baseline to a specified value for later use in q_hat.

Parameters:d (dict) – Dictionary containing data to insert into data-weighting R matrix cache. Keys are tuples, following the same format as the input to self.R().
set_iC(d)[source]

Set the cached inverse covariance matrix for a given dataset and baseline to a specified value. For now, you should already have applied weights to this matrix.

Parameters:d (dict) – Dictionary containing data to insert into inverse covariance matrix cache. Keys are tuples, following the same format as the input to self.iC().
set_spw(spw_range)[source]

Set the spectral window range.

Parameters:spw_range (tuple, contains start and end of spw in channel indices) – used to slice the frequency array
set_taper(taper)[source]

Set data tapering type.

Parameters:taper (str) – Type of data tapering. See aipy.dsp.gen_window for options.
set_weighting(data_weighting)[source]

Set data weighting type.

Parameters:data_weighting (str) – Type of data weightings. Options=[‘identity’, ‘iC’]
trim_dset_lsts(lst_tol=6)[source]

Assuming all datasets in self.dsets are locked to the same LST grid (but each may have a constant offset), trim LSTs from each dset that aren’t found in all other dsets (within some decimal tolerance specified by lst_tol).

Warning: this edits the data in dsets in-place, and is not reversible.

Parameters:lst_tol (float) – Decimal tolerance [radians] for comparing float-valued LST bins.
units(little_h=True)[source]

Return the units of the power spectrum. These are inferred from the units reported by the input visibilities (UVData objects).

Parameters:little_h (boolean, optional) – Whether to have cosmological length units be h^-1 Mpc or Mpc Default: h^-1 Mpc
Returns:pspec_units – Units of the power spectrum that is returned by pspec().
Return type:str
validate_datasets(verbose=True)[source]

Validate stored datasets and weights to make sure they are consistent with one another (e.g. have the same shape, baselines etc.).

validate_pol(dsets, pol_pair)[source]

Validate polarization and returns the index of the datasets so that the polarization pair is consistent with the UVData objects.

Parameters:
  • dsets (length-2 list or length-2 tuple of integers or str) – Contains indices of self.dsets to use in forming power spectra, where the first index is for the Left-Hand dataset and second index is used for the Right-Hand dataset (see above).
  • pol_pair (length-2 tuple of integers or strings) – Contains polarization pair which will be used in estiamting the power spectrum e,g (-5, -5) or (‘xy’, ‘xy’). Only auto-polarization pairs are implemented for the time being.
Returns:

valid – True if the UVData objects polarizations are consistent with the pol_pair (user specified polarizations) else False.

Return type:

boolean

w(key)[source]

Get weights for a given dataset and baseline, as specified in a standard key format.

Parameters:key (tuple) – Tuple containing dataset ID and baseline index. The first element of the tuple is the dataset index, and the subsequent elements are the baseline ID.
Returns:w – Array of weights for the requested UVData dataset and baseline.
Return type:array_like
x(key)[source]

Get data for a given dataset and baseline, as specified in a standard key format.

Parameters:key (tuple) – Tuple containing dataset ID and baseline index. The first element of the tuple is the dataset index (or label), and the subsequent elements are the baseline ID.
Returns:x – Array of data from the requested UVData dataset and baseline.
Return type:array_like