UVData Class

UVData is the main user class. It provides import and export functionality to all supported file formats (UVFITS, Miriad, FHD) and can be interacted with directly.

class pyuvdata.UVData[source]

A class for defining a radio interferometer dataset.

Currently supported file types: uvfits, miriad, fhd. Provides phasing functions.

UVParameter objects

For full list see UVData Parameters (http://pyuvdata.readthedocs.io/en/latest/uvdata_parameters.html). Some are always required, some are required for certain phase_types and others are always optional.

antnums_to_baseline(ant1, ant2, attempt256=False)[source]

Get the baseline number corresponding to two given antenna numbers.

Parameters:
  • ant1 – first antenna number (integer)
  • ant2 – second antenna number (integer)
  • attempt256 – Option to try to use the older 256 standard used in many uvfits files (will use 2048 standard if there are more than 256 antennas). Default is False.
Returns:

integer baseline number corresponding to the two antenna numbers.

antpair2ind(ant1, ant2)[source]

Get blt indices for given (ordered) antenna pair.

antpairpol_iter(squeeze='default')[source]

Generates numpy arrays of data for each antpair, pol combination.

Parameters:squeeze – ‘default’: squeeze pol and spw dimensions if possible ‘none’: no squeezing of resulting numpy array ‘full’: squeeze all length 1 dimensions
Returns (for each iteration):
key: tuple with antenna1, antenna2, and polarization string data: Numpy array with data which is the result of self[key]
baseline_to_antnums(baseline)[source]

Get the antenna numbers corresponding to a given baseline number.

Parameters:baseline – integer baseline number
Returns:tuple with the two antenna numbers corresponding to the baseline.
check(check_extra=True, run_check_acceptability=True)[source]

Add some extra checks on top of checks on UVBase class.

Check that required parameters exist. Check that parameters have appropriate shapes and optionally that the values are acceptable.

Parameters:
  • check_extra – If true, check all parameters, otherwise only check required parameters.
  • run_check_acceptability – Option to check if values in parameters are acceptable. Default is True.
get_ENU_antpos(center=True, pick_data_ants=False)[source]

Returns antenna positions in ENU (topocentric) coordinates in units of meters.

center : bool, if True: subtract median of array position from antpos pick_data_ants : bool, if True: return only antennas found in data

antpos : ndarray, antenna positions in TOPO frame and units of meters, shape=(Nants, 3) ants : ndarray, antenna numbers matching ordering of antpos, shape=(Nants,)

get_antpairpols()[source]

Returns list of unique antpair + pol tuples (ant1, ant2, pol) in data.

get_antpairs()[source]

Returns list of unique antpair tuples (ant1, ant2) in data.

get_ants()[source]

Returns numpy array of unique antennas in the data.

get_baseline_nums()[source]

Returns numpy array of unique baseline numbers in data.

get_data(*args, **kwargs)[source]

Function for quick access to numpy array with data corresponding to a baseline and/or polarization. Returns a read-only view if possible, otherwise a copy.

Parameters:
  • *args – parameters or tuple of parameters defining the key to identify desired data. See _key2inds for formatting.
  • **kwargs

    Keyword arguments: squeeze: ‘default’: squeeze pol and spw dimensions if possible

    ’none’: no squeezing of resulting numpy array ‘full’: squeeze all length 1 dimensions
    force_copy: Option to explicitly make a copy of the data.
    Default is False.
Returns:

Numpy array of data corresponding to key. If data exists conjugate to requested antenna pair, it will be conjugated before returning.

get_feedpols()[source]

Returns list of unique antenna feed polarizations (e.g. [‘X’, ‘Y’]) in data. NOTE: Will return ValueError if any pseudo-Stokes visibilities are present.

get_flags(*args, **kwargs)[source]

Function for quick access to numpy array with flags corresponding to a baseline and/or polarization. Returns a read-only view if possible, otherwise a copy.

Parameters:
  • *args – parameters or tuple of parameters defining the key to identify desired data. See _key2inds for formatting.
  • **kwargs

    Keyword arguments: squeeze: ‘default’: squeeze pol and spw dimensions if possible

    ’none’: no squeezing of resulting numpy array ‘full’: squeeze all length 1 dimensions
    force_copy: Option to explicitly make a copy of the data.
    Default is False.
Returns:

Numpy array of flags corresponding to key.

get_nsamples(*args, **kwargs)[source]

Function for quick access to numpy array with nsamples corresponding to a baseline and/or polarization. Returns a read-only view if possible, otherwise a copy.

Parameters:
  • *args – parameters or tuple of parameters defining the key to identify desired data. See _key2inds for formatting.
  • **kwargs

    Keyword arguments: squeeze: ‘default’: squeeze pol and spw dimensions if possible

    ’none’: no squeezing of resulting numpy array ‘full’: squeeze all length 1 dimensions
    force_copy: Option to explicitly make a copy of the data.
    Default is False.
Returns:

Numpy array of nsamples corresponding to key.

get_pols()[source]

Returns list of pols in string format.

get_times(*args)[source]

Find the time_array entries for a given antpair or baseline number. Meant to be used in conjunction with get_data function.

Parameters:*args – parameters or tuple of parameters defining the key to identify desired data. See _key2inds for formatting.
Returns:Numpy array of times corresonding to key.
juldate2ephem(num)[source]

Convert Julian date to ephem date, measured from noon, Dec. 31, 1899.

Parameters:num – Julian date
Returns:ephem date, measured from noon, Dec. 31, 1899.
known_telescopes()[source]

Retun a list of telescopes known to pyuvdata.

This is just a shortcut to uvdata.telescopes.known_telescopes()

order_pols(order='AIPS')[source]

Arranges polarizations in orders corresponding to either AIPS or CASA convention. CASA stokes types are ordered with cross-pols in between (e.g. XX,XY,YX,YY) while AIPS orders pols with auto-pols followed by cross-pols (e.g. XX,YY,XY,YX) Args: order: string, either ‘CASA’ or ‘AIPS’. Default=’AIPS’

parse_ants(ant_str, print_toggle=False)[source]

Generates two lists of antenna pair tuples and polarization indices based on parsing of the string ant_str. If no valid polarizations (pseudo-Stokes params, or combinations of [lr] or [xy]) or antenna numbers are found in ant_str, ant_pairs_nums and polarizations are returned as None.

Parameters:
  • ant_str – String containing antenna information to pass to select function. Can be ‘all’, ‘auto’, cross, or combinations of antenna numbers and polarization indicators l and r or x and y. Minus signs can also be used in front of an antenna number or baseline to exclude it from being output in ant_pairs_nums. If the antenna number attached with a minus sign is present in the outputted list ant_pairs_nums, it will be removed from the list. If ant_str passed with a minus sign as the first character, ‘all,’ will be appended to the beginning of the string. See the tutorial for examples of valid strings and their behavior.
  • print_toggle – Boolean for printing parsed baselines for a visual user check.
Output:
ant_pairs_nums: List of tuples containing the parsed pairs of
antenna numbers. If ‘all’ or pseudo-Stokes parameters are passed as ant_str, returned as None.
polarizations: List of desired polarizations. If no polarizations
found in ant_str then returned as None.
phase(ra, dec, epoch)[source]

Phase a drift scan dataset to a single ra/dec at a particular epoch.

Will not phase already phased data.

Parameters:
  • ra – The ra to phase to in radians.
  • dec – The dec to phase to in radians.
  • epoch – The epoch to use for phasing. Should be an ephem date, measured from noon Dec. 31, 1899.
phase_to_time(time)[source]

Phase a drift scan dataset to the ra/dec of zenith at a particular time.

Parameters:time – The time to phase to.
read_fhd(filelist, use_model=False, run_check=True, check_extra=True, run_check_acceptability=True)[source]

Read in data from a list of FHD files.

Parameters:
  • filelist – The list of FHD save files to read from. Must include at least one polarization file, a params file and a flag file. Can also be a list of lists to read multiple data sets.
  • use_model – Option to read in the model visibilities rather than the dirty visibilities. Default is False.
  • run_check – Option to check for the existence and proper shapes of parameters after reading in the file. Default is True.
  • check_extra – Option to check optional parameters as well as required ones. Default is True.
  • run_check_acceptability – Option to check acceptable range of the values of parameters after reading in the file. Default is True.
read_miriad(filepath, correct_lat_lon=True, run_check=True, check_extra=True, run_check_acceptability=True, phase_type=None, antenna_nums=None, ant_str=None, ant_pairs_nums=None, polarizations=None, time_range=None)[source]

Read in data from a miriad file.

Parameters:
  • filepath – The miriad file directory or list of directories to read from.
  • run_check – Option to check for the existence and proper shapes of parameters after reading in the file. Default is True.
  • check_extra – Option to check optional parameters as well as required ones. Default is True.
  • run_check_acceptability – Option to check acceptable range of the values of parameters after reading in the file. Default is True.
  • antenna_nums – The antennas numbers to read into the object.
  • ant_pairs_nums – A list of antenna number tuples (e.g. [(0,1), (3,2)]) specifying baselines to read into the object. Ordering of the numbers within the tuple does not matter. A single antenna iterable e.g. (1,) is interpreted as all visibilities with that antenna.
  • ant_str – A string containing information about what kinds of visibility data to read-in. Can be ‘auto’, ‘cross’, ‘all’. Cannot provide ant_str if antenna_nums and/or ant_pairs_nums is not None.
  • polarizations – List of polarization integers or strings to read-in. Ex: [‘xx’, ‘yy’, …]
  • time_range – len-2 list containing min and max range of times (Julian Date) to read-in. Ex: [2458115.20, 2458115.40]
read_ms(filepath, run_check=True, check_extra=True, run_check_acceptability=True, data_column='DATA', pol_order='AIPS')[source]

Read in data from a measurement set

Parameters:
  • filepath – The measurement set file directory or list of directories to read from.
  • run_check – Option to check for the existence and proper shapes of parameters after reading in the file. Default is True.
  • check_extra – Option to check optional parameters as well as required ones. Default is True.
  • run_check_acceptability – Option to check the values of parameters after reading in the file. Default is True.
  • data_column – name of CASA data column to read into data_array. ‘DATA’, ‘MODEL’, or ‘CORRECTED_DATA’
  • pol_order – specify whether you want polarizations ordered by ‘CASA’ or ‘AIPS’ conventions.
read_uvfits(filename, antenna_nums=None, antenna_names=None, ant_str=None, ant_pairs_nums=None, frequencies=None, freq_chans=None, times=None, polarizations=None, blt_inds=None, read_data=True, read_metadata=True, run_check=True, check_extra=True, run_check_acceptability=True)[source]

Read in header, metadata and data from uvfits file(s).

Parameters:
  • filename – The uvfits file or list of files to read from.
  • antenna_nums – The antennas numbers to include when reading data into the object (antenna positions and names for the excluded antennas will be retained). This cannot be provided if antenna_names is also provided. Ignored if read_data is False.
  • antenna_names – The antennas names to include when reading data into the object (antenna positions and names for the excluded antennas will be retained). This cannot be provided if antenna_nums is also provided. Ignored if read_data is False.
  • ant_pairs_nums – A list of antenna number tuples (e.g. [(0,1), (3,2)]) specifying baselines to include when reading data into the object. Ordering of the numbers within the tuple does not matter. Ignored if read_data is False.
  • ant_str – A string containing information about what antenna numbers and polarizations to include when reading data into the object. Can be ‘auto’, ‘cross’, ‘all’, or combinations of antenna numbers and polarizations (e.g. ‘1’, ‘1_2’, ‘1x_2y’). See tutorial for more examples of valid strings and the behavior of different forms for ant_str. If ‘1x_2y,2y_3y’ is passed, both polarizations ‘xy’ and ‘yy’ will be kept for both baselines (1,2) and (2,3) to return a valid pyuvdata object. An ant_str cannot be passed in addition to any of the above antenna args or the polarizations arg. Ignored if read_data is False.
  • frequencies – The frequencies to include when reading data into the object. Ignored if read_data is False.
  • freq_chans – The frequency channel numbers to include when reading data into the object. Ignored if read_data is False.
  • times – The times to include when reading data into the object. Ignored if read_data is False.
  • polarizations – The polarizations to include when reading data into the object. Ignored if read_data is False.
  • blt_inds – The baseline-time indices to include when reading data into the object. This is not commonly used. Ignored if read_data is False.
  • read_data – Read in the visibility and flag data. If set to false, only the basic header info and metadata (if read_metadata is True) will be read in. Results in an incompletely defined object (check will not pass). Default True.
  • read_metadata – Read in metadata (times, baselines, uvws) as well as basic header info. Only used if read_data is False (metadata will be read if data is read). If both read_data and read_metadata are false, only basic header info is read in. Default True.
  • run_check – Option to check for the existence and proper shapes of parameters after reading in the file. Default is True. Ignored if read_data is False.
  • check_extra – Option to check optional parameters as well as required ones. Default is True. Ignored if read_data is False.
  • run_check_acceptability – Option to check acceptable range of the values of parameters after reading in the file. Default is True. Ignored if read_data is False.
read_uvfits_data(filename, antenna_nums=None, antenna_names=None, ant_str=None, ant_pairs_nums=None, frequencies=None, freq_chans=None, times=None, polarizations=None, blt_inds=None, run_check=True, check_extra=True, run_check_acceptability=True)[source]

Read in data but not header info from a uvfits file (useful for an object that already has the associated header info).

Parameters:
  • filename – The uvfits file
  • antenna_nums – The antennas numbers to include when reading data into the object (antenna positions and names for the excluded antennas will be retained). This cannot be provided if antenna_names is also provided.
  • antenna_names – The antennas names to include when reading data into the object (antenna positions and names for the excluded antennas will be retained). This cannot be provided if antenna_nums is also provided.
  • ant_pairs_nums – A list of antenna number tuples (e.g. [(0,1), (3,2)]) specifying baselines to include when reading data into the object. Ordering of the numbers within the tuple does not matter.
  • ant_str – A string containing information about what antenna numbers and polarizations to include when reading data into the object. Can be ‘auto’, ‘cross’, ‘all’, or combinations of antenna numbers and polarizations (e.g. ‘1’, ‘1_2’, ‘1x_2y’). See tutorial for more examples of valid strings and the behavior of different forms for ant_str. If ‘1x_2y,2y_3y’ is passed, both polarizations ‘xy’ and ‘yy’ will be kept for both baselines (1,2) and (2,3) to return a valid pyuvdata object. An ant_str cannot be passed in addition to any of the above antenna args or the polarizations arg.
  • frequencies – The frequencies to include when reading data into the object.
  • freq_chans – The frequency channel numbers to include when reading data into the object.
  • times – The times to include when reading data into the object.
  • polarizations – The polarizations to include when reading data into the object.
  • blt_inds – The baseline-time indices to include when reading data into the object. This is not commonly used.
  • run_check – Option to check for the existence and proper shapes of parameters after reading in the file. Default is True.
  • check_extra – Option to check optional parameters as well as required ones. Default is True.
  • run_check_acceptability – Option to check acceptable range of the values of parameters after reading in the file. Default is True.
read_uvfits_metadata(filename)[source]

Read in metadata (random parameter info) but not data from a uvfits file (useful for an object that already has the associated header info and full visibility data isn’t needed).

Parameters:filename – The uvfits file to read from.
read_uvh5(filename, run_check=True, check_extra=True, run_check_acceptability=True)[source]

Read a UVH5 file.

Parameters:
  • filename – The UVH5 file to read.
  • run_check – Option to check for the existence and proper shapes of parameters after reading in the file. Default is True.
  • check_extra – Option to check optional parameters as well as required ones. Default is True.
  • run_check_acceptability – Option to check acceptable range of the values of parameters after reading in the file. Default is True.
Returns:

None

reorder_pols(order=None, run_check=True, check_extra=True, run_check_acceptability=True)[source]

Rearrange polarizations in the event they are not uvfits compatible.

Parameters:
  • order – Provide the order which to shuffle the data. Default will sort by absolute value of pol values.
  • run_check – Option to check for the existence and proper shapes of parameters after reordering. Default is True.
  • check_extra – Option to check optional parameters as well as required ones. Default is True.
  • run_check_acceptability – Option to check acceptable range of the values of parameters after reordering. Default is True.
select(antenna_nums=None, antenna_names=None, ant_str=None, ant_pairs_nums=None, frequencies=None, freq_chans=None, times=None, polarizations=None, blt_inds=None, run_check=True, check_extra=True, run_check_acceptability=True, inplace=True)[source]

Select specific antennas, antenna pairs, frequencies, times and polarizations to keep in the object while discarding others.

Also supports selecting specific baseline-time indices to keep while discarding others, but this is not commonly used. The history attribute on the object will be updated to identify the operations performed.

Parameters:
  • antenna_nums – The antennas numbers to keep in the object (antenna positions and names for the removed antennas will be retained). This cannot be provided if antenna_names is also provided.
  • antenna_names – The antennas names to keep in the object (antenna positions and names for the removed antennas will be retained). This cannot be provided if antenna_nums is also provided.
  • ant_pairs_nums – A list of antenna number tuples (e.g. [(0,1), (3,2)]) specifying baselines to keep in the object. Ordering of the numbers within the tuple does not matter.
  • ant_str – A string containing information about what antenna numbers and polarizations to keep in the object. Can be ‘auto’, ‘cross’, ‘all’, or combinations of antenna numbers and polarizations (e.g. ‘1’, ‘1_2’, ‘1x_2y’). See tutorial for more examples of valid strings and the behavior of different forms for ant_str. If ‘1x_2y,2y_3y’ is passed, both polarizations ‘xy’ and ‘yy’ will be kept for both baselines (1,2) and (2,3) to return a valid pyuvdata object. An ant_str cannot be passed in addition to any of the above antenna args or the polarizations arg. If this occurs, select raises a ValueError.
  • frequencies – The frequencies to keep in the object.
  • freq_chans – The frequency channel numbers to keep in the object.
  • times – The times to keep in the object.
  • polarizations – The polarizations to keep in the object.
  • blt_inds – The baseline-time indices to keep in the object. This is not commonly used.
  • run_check – Option to check for the existence and proper shapes of parameters after downselecting data on this object. Default is True.
  • check_extra – Option to check optional parameters as well as required ones. Default is True.
  • run_check_acceptability – Option to check acceptable range of the values of parameters after downselecting data on this object. Default is True.
  • inplace – Option to perform the select directly on self (True, default) or return a new UVData object, which is a subselection of self (False)
set_drift()[source]

Set phase_type to ‘drift’ and adjust required parameters.

set_lsts_from_time_array()[source]

Set the lst_array based from the time_array.

set_phased()[source]

Set phase_type to ‘phased’ and adjust required parameters.

set_telescope_params(overwrite=False)[source]

Set telescope related parameters.

If the telescope_name is in the known_telescopes, set any missing telescope-associated parameters (e.g. telescope location) to the value for the known telescope.

Parameters:overwrite – Option to overwrite existing telescope-associated parameters with the values from the known telescope. Default is False.
set_unknown_phase_type()[source]

Set phase_type to ‘unknown’ and adjust required parameters.

unphase_to_drift()[source]

Convert from a phased dataset to a drift dataset.

write_miriad(filepath, run_check=True, check_extra=True, run_check_acceptability=True, clobber=False, no_antnums=False)[source]

Write the data to a miriad file.

Parameters:
  • filename – The miriad file directory to write to.
  • run_check – Option to check for the existence and proper shapes of parameters before writing the file. Default is True.
  • check_extra – Option to check optional parameters as well as required ones. Default is True.
  • run_check_acceptability – Option to check acceptable range of the values of parameters before writing the file. Default is True.
  • clobber – Option to overwrite the filename if the file already exists. Default is False.
  • no_antnums – Option to not write the antnums variable to the file. Should only be used for testing purposes.
write_uvfits(filename, spoof_nonessential=False, force_phase=False, run_check=True, check_extra=True, run_check_acceptability=True)[source]

Write the data to a uvfits file.

Parameters:
  • filename – The uvfits file to write to.
  • spoof_nonessential – Option to spoof the values of optional UVParameters that are not set but are required for uvfits files. Default is False.
  • force_phase – Option to automatically phase drift scan data to zenith of the first timestamp. Default is False.
  • run_check – Option to check for the existence and proper shapes of parameters before writing the file. Default is True.
  • check_extra – Option to check optional parameters as well as required ones. Default is True.
  • run_check_acceptability – Option to check acceptable range of the values of parameters before writing the file. Default is True.
write_uvh5(filename, run_check=True, check_extra=True, run_check_acceptability=True, clobber=False)[source]

Write a UVData object to a UVH5 file.

Parameters:
  • filename – The UVH5 file to write to.
  • run_check – Option to check for the existence and proper shapes of parameters before writing the file. Default is True.
  • check_extra – Option to check optional parameters as well as required ones. Default is True.
  • run_check_acceptability – Option to check acceptable range of the values of parameters before writing the file. Default is True.
  • clobber – Option to overwrite the file if it already exists. Default is False.
Returns:

None