Imaging

Data Representation

Classes to hold visibility data that are used internally with LSL.

class lsl.imaging.data.PolarizationDataSet(polarization, data=None, weight=None, mask=None)

Class to hold the visibility data/weight/mask for a single polarization.

append(data, weight=None, mask=None)

Append a new visibility to the object.

copy()

Return a copy of the object.

extend(data, weight=None, mask=None)

Extend the object with a collection of visibilites.

property nbaseline

The number of baselines contained in the data.

sort(order)

Apply the specified ordering on the object.

subset(selection)

Return a subset of the data.

class lsl.imaging.data.VisibilityData(data=None)

Class to hold multiple integrations of visibility data.

append(value)

Append a new integration stored as a VisibilityDataSet to the object.

property baselines

Return a list baselines contained in the first VisbilityDataSet object.

copy(include_pols=True)

Return a copy of the object. Be default this includes copies of all of the associated PolarizationDataSet objects stored in the individual VisibilityDataSet sets. Setting ‘include_pols’ to False will not copy these objects.

extend(values)

Append a collection of new integration stored as VisibilityDataSet to the objects.

property freq

Return a list polarizations contained in the first VisbilityDataSet object.

get_antenna_subset(include=None, exclude=None, indicies=True)

Return a copy of the data containing baselines that meet the specified antenna selection criteria. The selection is governed by the ‘include’, ‘exclude’, and ‘indicies’ keywords. If ‘include’ is not None, only baselines containing antennas in the list are selected. If ‘exclude’ is not None, only baselines not containing antennas in the list are selected. If both ‘include’ and ‘exclude’ are provided, the ‘include’ list has priority. ‘indicies’ sets whether or not the selection criteria are provided as indicies of the antennas stored in the antennaarray attribute or are stand numbers.

get_uv_range(min_uv=0.0, max_uv=inf)

Return a copy of the data containing only baselines with the (u,v) distances allowed by the ‘min_uv’ and ‘max_uv’ cuts.

property jds

Return a list of JDs for all VisibilityDataSets contained.

property mjds

Return a list of MJDs for all VisibilityDataSets contained.

property nbaseline

The number of baselines contained in the first VisbilityDataSet object.

property nchan

The number of frequency channels contained in the first VisbilityDataSet object.

property npol

The number of polarizations contained in the first VisbilityDataSet object.

property pols

Return a list polarizations contained in the first VisbilityDataSet object.

pop(index=- 1)

Pop and return the VisibilityDataSet specified by the provided index.

rephase(new_phase_center, ignore_errors=False)

Shift the phase center for all of the integrations stored.

sort()

Sort the VisibilityDataSet objects contained by the JD of the integrations.

class lsl.imaging.data.VisibilityDataSet(jd, freq, baselines, uvw, antennaarray=None, phase_center=None, **kwds)

Class to wrap a single integration of visibility data.

append(value)

Append a new PolarizationDataSet object. If the polarization already exists it is replaced.

copy(include_pols=True)

Return a copy of the object. Be default this includes copies of all of the associated PolarizationDataSet objects. Setting ‘include_pols’ to False will not copy these objects.

get_antenna_subset(include=None, exclude=None, indicies=True)

Return a copy of the data containing baselines that meet the specified antenna selection criteria. The selection is governed by the ‘include’, ‘exclude’, and ‘indicies’ keywords. If ‘include’ is not None, only baselines containing antennas in the list are selected. If ‘exclude’ is not None, only baselines not containing antennas in the list are selected. If both ‘include’ and ‘exclude’ are provided, the ‘include’ list has priority. ‘indicies’ sets whether or not the selection criteria are provided as indicies of the antennas stored in the antennaarray attribute or are stand numbers.

get_uv_range(min_uv=0.0, max_uv=inf)

Return a copy of the data containing only baselines with the (u,v) distances allowed by the ‘min_uv’ and ‘max_uv’ cuts.

property mjd

The MJD that the data correspond to.

property nbaseline

The number of baselines in the data.

property nchan

The number of frequency channels in the data.

property npol

The number of polarizations stored.

rephase(new_phase_center)

Shift the phase center of the data to a new phase center.

sort(order=None)

Sort the stored data using the order provided in the ‘order’ keyword. If not ordering is provided, the data are sorted by baseline pair.

Imaging

Module to support imaging correlated data. This module provides utilities to read FITS IDI files into lsl.imaging.data.VisibilityDataSet or lsl.imaging.data.VisibilityData object (as described in lsl.imaging.data) and build AIPY ImgW instances from the data.

New in version 0.5.0.

Changed in version 1.0.0: Added support for UVFITS files and CASA measurement sets

Changed in version 1.0.1: Added the plot_gridded_image() function

Changed in version 1.1.0: Added the get_image_radec() and get_image_azalt() functions to complement plot_gridded_image() and make it easier to work with phase centers that are not at zenith. Added in the ImgWPlus class to add support for imaging weighting and tapering.

Changed in version 1.1.2: Added support for CASA MeasurementSets that are stored as a tarball

Changed in version 1.1.4: Fixed a conjugation problem in the visibilities read from a FITS-IDI file

lsl.imaging.utils.CorrelatedData(filename, verbose=False)

Read in and work with FITS IDI and UVFITS files. Returns either a CorrelateDataIDI or CorrelatedDataUV instance.

class lsl.imaging.utils.CorrelatedDataIDI(filename)

Class to make accessing information about a FITS IDI easy. This wraps all of the “messy” machinery needed to extract both the metadata and data from the file and return them as common LSL objects.

This class has three main attributes to interact with:
  • get_antennaarray - Return a lsl.sim.vim.AntennaArray instance

    that represents the array where the data was obtained. This is useful for simulation proposes and computing source positions.

  • get_observer - Return a ephem.Observer instance representing the array

  • get_data_set - Return a lsl.imaging.data.VisibilityDataSet or

    lsl.imaging.data.VisibilityData object for all baselines for a given set of observations

The class also includes a variety of useful metadata attributes:
  • pols - Numpy array of polarization product codes

  • freq - Numpy array of frequency channels in Hz

  • station - LSL lsl.common.stations.LWAStation instance for the

    array

  • date_obs - Datetime object for the reference date of the FIT IDI file

  • antennas - List of lsl.common.stations.Antenna instances

Note

The CorrelatedData.antennas attribute should be used over CorrelatedData.station.antennas since the mapping in the FITS IDI file may not be the same as the digitizer order.

get_data_set(sets, include_auto=False, sort=True, min_uv=0, max_uv=inf)

Return a lsl.imaging.data.VisibilityDataSet or lsl.imaging.data.VisibilityData object for all baselines for a given set of observations for the specified data set. By default this excludes the autocorrelations. To include autocorrelations set the value of ‘include_auto’ to True. Setting the ‘sort’ keyword to False will disable the baseline sorting. Optionally, baselines with lengths between min_uv and max_uv can only be returned.

Note

min_uv and max_uv should be specified in lambda

Changed in version 1.1.0: ‘set’ can now be either an integer or a list to pull back multiple integrations.

class lsl.imaging.utils.CorrelatedDataMS(filename)

Class to make accessing information about a MS easy. This wraps all of the “messy” machinery needed to extract both the metadata and data from the file and return them as common LSL objects.

class lsl.imaging.utils.CorrelatedDataUV(filename)

Class to make accessing information about a UVFITS file easy. This wraps all of the “messy” machinery needed to extract both the metadata and data from the file and return them as common LSL objects.

This class has three main attributes to interact with:
  • get_antennaarray - Return a lsl.sim.vim.AntennaArray instance

    that represents the array where the data was obtained. This is useful for simulation proposes and computing source positions.

  • get_observer - Return a ephem.Observer instance representing the array

  • get_data_set - Return a lsl.imaging.data.VisibilityDataSet or

    lsl.imaging.data.VisibilityData object of all baselines for a given set of observations

The class also includes a variety of useful metadata attributes:
  • pols - Numpy array of polarization product codes

  • freq - Numpy array of frequency channels in Hz

  • station - LSL lsl.common.stations.LWAStation instance for the

    array

  • date_obs - Datetime object for the reference date of the FIT IDI file

  • antennas - List of lsl.common.stations.Antenna instances

Note

The CorrelatedDataUV.antennas attribute should be used over CorrelatedDataUV.station.antennas since the mapping in the UVFITS file may not be the same as the digitizer order.

get_data_set(sets, include_auto=False, sort=True, min_uv=0, max_uv=inf)

Return a lsl.imaging.data.VisibilityDataSet or lsl.imaging.data.VisibilityData object for the specified data set. By default this excludes the autocorrelations. To include autocorrelations set the value of ‘include_auto’ to True. Setting the ‘sort’ keyword to False will disable the baseline sorting. Optionally, baselines with lengths between min_uv and max_uv can only be returned.

Note

min_uv and max_uv should be specified in lambda

Changed in version 1.1.0: ‘set’ can now be either an integer or a list to pull back multiple integrations.

class lsl.imaging.utils.ImgWPlus(size=100, res=1, wres=0.5, mf_order=0)

Sub-class of the aipy.img.ImgW class that adds support for different visibility weighting scheme and uv plane tapering. This class also adds in a couple of additional methods that help determine the size of the field of view and the pixels near the phase center.

bm_image(center=0, 0, term=None, weighting='natural', local_fraction=0.5, robust=0.0, taper=0.0, 0.0)

Return the inverse FFT of the sample weightings (for all mf_order terms, or the specified term if supplied), with the 0,0 point moved to ‘center’. In the images return north is up and east is to the left.

There are a few keywords that control how the image is formed. There are:

  • weighting - The weighting scheme (‘natural’, ‘uniform’, or

    ‘briggs’) used on the data;

  • local_fraction - The fraction of the uv grid that is consider

    “local” for the ‘uniform’ and ‘briggs’ methods;

  • robust - The value for the weighting robustness under the

    ‘briggs’ method; and

  • taper - The size of u and v Gaussian tapers at the 30% level.

property field_of_view

Return the approximate size of the field of view in radians. The field of view calculate is based off the maximum and minimum values of L found for the inverted uv matrix.

image(center=0, 0, weighting='natural', local_fraction=0.5, robust=0.0, taper=0.0, 0.0)

Return the inverse FFT of the UV matrix, with the 0,0 point moved to ‘center’. In the images return north is up and east is to the left.

There are a few keywords that control how the image is formed. There are:

  • weighting - The weighting scheme (‘natural’, ‘uniform’, or

    ‘briggs’) used on the data;

  • local_fraction - The fraction of the uv grid that is consider

    “local” for the ‘uniform’ and ‘briggs’ methods;

  • robust - The value for the weighting robustness under the

    ‘briggs’ method; and

  • taper - The size of u and v Gaussian tapers at the 30% level.

property pixel_size

Return the approximate size of pixels at the phase center in radians. The pixel size is averaged over the four pixels that neighboor the phase center.

put(uvw, data, wgts=None, invker2=None, verbose=True)

Same as Img.put, only now the w component is projected to the w=0 plane before applying the data to the UV matrix.

lsl.imaging.utils.build_gridded_image(data_set, size=80, res=0.5, wres=0.1, pol='XX', chan=None, im=None, verbose=True)

Given a lsl.imaging.data.VisibilityDataSet object, build an aipy.img.ImgW object of gridded uv data which can be used for imaging. The ImgW object itself is returned by this function to make it more versatile.

lsl.imaging.utils.get_image_azalt(gimg, aa, phase_center='z', shifted=True)

Given a gridded image generated by the build_gridded_image() function and an AntennaArray instance, return a two-element tuple containing the azimuth and elevation (altitude), both in radians, for each pixel in the image.

The ‘phase_center’ keyword controls what the phase center of the image is and defaults to zenith.

lsl.imaging.utils.get_image_radec(gimg, aa, phase_center='z', shifted=True)

Given a gridded image generated by the build_gridded_image() function and an AntennaArray instance, return a two-element tuple containing the RA and dec. values (in radians) for each pixel in the image.

The ‘phase_center’ keyword controls what the phase center of the image is and defaults to zenith.

lsl.imaging.utils.plot_gridded_image(ax, gimg, shifted=True, origin='lower', interpolation='nearest', **kwargs)

Given a blank matplotlib axes instance and a gridded image generated by the build_gridded_image() function, plot the image on the axes and setup the basic coordinate system. This function returns the matplotlib object added to the plot

Changed in version 1.2.1: Changed the function to return the matplotlib object plotted so that colorbars can be added

Changed in version 1.1.0: Added a ‘shifted’ keyword to control whether or not the image is centered or not.

New in version 1.0.1.

Plot Overlays

Module that provides a variety of overlays for all-sky images. These overlays include:

  • the locations and names of sources,

  • the horizon,

  • a graticle showing lines of constant RA and dec., and

  • a graticle showing lines of constant azimuth and elevation.

All of the functions in this module accept a matplotlib axes instances that is used for plotting.

New in version 1.0.1.

Changed in version 1.1.0: Added support for overlaying on images with non-zenith phase centers

lsl.imaging.overlay.graticule_azalt(ax, antennaarray, phase_center='z', label=True, color='white')

For a matplotlib axis instance showing an image of the sky, plot lines of constant azimuth and elevation. Elevations are spaced at 20 degree intervals and azimuths are spaced at 45 degree intervals

lsl.imaging.overlay.graticule_radec(ax, antennaarray, phase_center='z', label=True, color='white')

For a matplotlib axis instance showing an image of the sky, plot lines of constant declinate and RA. Declinations are spaced at 20 degree intervals and RAs are spaced at 2 hour intervals.

lsl.imaging.overlay.horizon(ax, antennaarray, phase_center='z', color='white')

For a matplotlib axis instance showing an image of the sky, plot the horizon.

Changed in version 1.1.0: Added a new argument for the AntennaArray instance to provide a uniform default call for all functions.

lsl.imaging.overlay.sources(ax, antennaarray, srcs, phase_center='z', label=True, marker='x', color='white')

For a matplotlib axis instance showing an image of the sky, plot the locations of the srcs given in the ‘srcs’ dictionary.

Basic Image Analysis

Module for analyzing images. Currently, this module supports:
  • estimating the position-dependent background level in the image

  • finding point sources

New in version 1.1.0.

lsl.imaging.analysis.estimate_background(image, window=32)

Given a 2-D image, estimate and return the background a la SExtractor.

This works by:
  1. Dividing the image into a number of half-overlapped tiles that are ‘window’ by ‘window’ pixels in size.

  2. Computing the mean and standard deviation within the tile both with and without 3*sigma clipping applied.

  3. Using the tile statistics, determine if the tile is empty, i.e., the standard deviation has changed by less than 20%, or full.

    1. If the tile is empty, use the clipped mean for the background.

    2. Otherwise, use 2.5*median - 1.5*mean.

  4. Once all of the tiles have been processes, median filter them to remove those that are anomalously high.

  5. Build a bicubic spline to interpolate the background back to the original image size.

  6. Evaluate the spline and return.

For more information on the SExtractor method, see:
lsl.imaging.analysis.find_point_sources(image, threshold=4.0, fwhm=1.0, sharp=[0.2, 1.0], round=[- 1.0, 1.0], background_size=16, verbose=True)

Given a 2-D image, find all of the point sources in it that meet the selection criteria provided via the keywords. These are:

  • threshold: detection threshold in counts about the background

  • fwhm: source full width half max. in pixel

  • sharp: two-element array that defines the lower and upper bounds

    on the source sharpness

  • round: two-element array that defines the lower and upper bounds

    on the source roundness

Background estimation and removal is handled by the estimate_background() function in this module that implements a SExtractor-like method. This can be disabled by setting the ‘background_size’ keyword to 0. For details see lsl.imaging.analysis.estimate_background().

The output of this function is a five-element tuple of 1-D NumPy arrays that store information for each source. The elements are:

  • x: intensity-weighted center - x coordinate

  • y: intensity-weighted center - y coordinate

  • flux: peak count value - input image units

  • sharpness: sharpness statistic

  • roundness: roundness statistic

This function is based on the FIND procedure from the AstroIDL User’s Library that was written by W. Landsman (Hughes STX).

For additional information about the original IDL routines, see: http://idlastro.gsfc.nasa.gov/contents.html#C2

Self-Calibration

Simple self-calibration module for correlated TBW and TBN data. The supported self-calibration methods are:

  • phase-only

  • amplitude and phase

  • delay-only

  • amplitude and delay

  • delay/phase offset

  • amplitude and delay/phase offset

..versionchanged:: 0.6.3

Reworked the module to make it more versatile

..versionadded:: 0.5.5

lsl.imaging.selfcal.delay_and_phase(aa, dataSet, simSet, chan, pol, ref_ant=0, max_iter=30, amplitude=False, amplitude_cutoff=1.001, delay_cutoff=0.2, phase_cutoff=0.01, verbose=True)

Function to apply a delay and phase offset (and, optionally, a amplitude) self-calibration to data stored in a readUVData dictionary and a model sky stored in a lsl.sim.vis.buildSimSky dictionary for the given polarization and channel(s).

Note

If the “amplitude” keyword is set to True, a four-element tuple of self-cal’d data, gain coefficients, delays in ns, and phase offsets is returned rather than the standard three-element tuple.

lsl.imaging.selfcal.delay_only(aa, dataSet, simSet, chan, pol, ref_ant=0, max_iter=30, amplitude=False, amplitude_cutoff=1.001, delay_cutoff=0.2, verbose=True)

Function to apply a delay-only (and, optionally, a amplitude) self- calibration to data stored in a readUVData dictionary and a model sky stored in a lsl.sim.vis.buildSimSky dictionary for the given polarization and channel(s).

Note

If the “amplitude” keyword is set to True, a three-element tuple of self-cal’d data, gain coefficients, and delays in ns is returned rather than the standard two-element tuple.

lsl.imaging.selfcal.phase_only(aa, dataSet, simSet, chan, pol, ref_ant=0, max_iter=30, amplitude=False, amplitude_cutoff=1.001, phase_cutoff=0.01, verbose=True)

Function to apply a phase-only (and, optionally, a amplitude) self- calibration to data stored in a readUVData dictionary and a model sky stored in a lsl.sim.vis.buildSimSky dictionary for the given polarization and channel(s).

Note

If the “amplitude” keyword is set to True, a three-element tuple of self-cal’d data, gain coefficients, and phase offsets is returned rather than the standard two-element tuple.

Deconvolution

Deconvolution support for images made with lsl.imaging.utils.build_gridded_image().

lsl.imaging.deconv.clean(aa, dataDict, aipyImg, input_image=None, size=80, res=0.5, wres=0.1, pol='XX', chan=None, gain=0.2, max_iter=150, sigma=3.0, verbose=True, plot=False)

Given a AIPY antenna array instance, a data dictionary, and an AIPY ImgW instance filled with data, return a deconvolved image. This function uses a CLEAN-like method that computes the array beam for each peak in the flux. Thus the CLEAN loop becomes:

  1. Find the peak flux in the residual image

  2. Compute the systems response to a point source at that location

  3. Remove the scaled porition of this beam from the residuals

  4. Go to 1.

CLEAN tuning parameters:
  • gain - CLEAN loop gain (default 0.2)

  • max_iter - Maximum number of iterations (default 150)

  • sigma - Threshold in sigma to stop cleaning (default 3.0)

lsl.imaging.deconv.clean_sources(aa, dataDict, aipyImg, srcs, input_image=None, size=80, res=0.5, wres=0.1, pol='XX', chan=None, gain=0.1, max_iter=150, sigma=2.0, verbose=True, plot=False)

Given a AIPY antenna array instance, a data dictionary, an AIPY ImgW instance filled with data, and a dictionary of sources, return the CLEAN components and the residuals map. This function uses a CLEAN-like method that computes the array beam for each peak in the flux. Thus the CLEAN loop becomes:

  1. Find the peak flux in the residual image

  2. Compute the systems response to a point source at that location

  3. Remove the scaled porition of this beam from the residuals

  4. Go to 1.

This function differs from clean() in that it only cleans localized regions around each source rather than the whole image. This is intended to help the mem() function along.

CLEAN tuning parameters:
  • gain - CLEAN loop gain (default 0.1)

  • max_iter - Maximum number of iterations (default 150)

  • sigma - Threshold in sigma to stop cleaning (default 2.0)

lsl.imaging.deconv.lsq(aa, dataDict, aipyImg, input_image=None, size=80, res=0.5, wres=0.1, pol='XX', chan=None, gain=0.05, max_iter=150, rtol=1e-09, verbose=True, plot=False)

Given a AIPY antenna array instance, a data dictionary, and an AIPY ImgW instance filled with data, return a deconvolved image. This function implements a least squares deconvolution.

Least squares tuning parameters:
  • gain - least squares loop gain (default 0.05)

  • max_iter - Maximum number of iteration (default 150)

  • rtol - Minimum change in the residual RMS between iterations

    (default 1e-9)