SKA-Low MCCS Calibration Package

High-level Setup and Calibration Functions.

ska_low_mccs_calibration.calibration.adjust_gsm(vis, modelvis, gal_data, sun_data)

Fit for a brightness scale difference between the galaxy and the sun.

The change may be polarised, but would need to peel the sun away for that.

Parameters:
  • vis (Dataset) – Calibrated Visibility xarray

  • modelvis (Dataset) – Model Visibility xarray to update

  • gal_data (ndarray[Any, dtype[Union[complex64, complex128]]]) – 2D NumPy ndarray containing galactic contribution to model vis

  • sun_data (ndarray[Any, dtype[Union[complex64, complex128]]]) – 2D NumPy ndarray containing solar contribution to model vis

Return type:

tuple[Dataset, float]

Returns:

(updated modelvis, adjustment factor)

ska_low_mccs_calibration.calibration.calibrate_mccs_visibility(vis, skymodel='gsm', masked_antennas=array([], dtype=int64), adjust_solar_model=True, jones_solve=False, nside=32, dtype=<class 'numpy.complex64'>, eep_path='.', eep_filebase='FEKO_AAVS3_vogel_256_elem_50ohm_', eep_rotation_deg=0.0)

Calibrate MCCS visibilities against a Global Sky Model.

Parameters:
  • vis (Dataset) – Visibility xarray to calibrate

  • skymodel (Union[Dataset, str]) – Either a model visibility xarray to calibrate against or the name of a sky model type to generate. The only supported sky model is “gsm”, which is the addition of any components above the horizon from the PyGDSM GlobalSkyModel16 Galactic sky model and a point-source solar model.

  • masked_antennas (ndarray, list) – list of antennas that have been flagged as bad and should be excluded. The default is to include all antennas.

  • adjust_solar_model (bool) – Adjust sky model after initial calibration to limit differences between the galactic and solar flux scales.

  • jones_solve (bool) – Call polarised solver after initial calibration.

  • nside (int) – Healpix nside at which to generate sky model. Default is 32, which corresponds to approximately 110 arcminutes. This is okay at low frequencies but may need to be increased at higher frequencies.

  • dtype (dtype) – The NumPy data-type to use for complex EEPs and sky models. The default is numpy.complex64.

  • eep_path (str) – Directory path to EEP image files. Default is the current dir.

  • eep_filebase (str) – Start of filenames, before frequency specification. Files are assumed to have a specific form, which for polarisation combination X/phi is: eep_path/eep_filebase<freq>MHz_Xpol_phi.npy. Default is filebase=”FEKO_AAVS3_vogel_256_elem_50ohm_”.

  • eep_rotation_deg (float) – Azimuthal rotation of station (and therefore the patterns) in degrees. Default is zero degress (no rotation).

Return type:

tuple[Dataset, Dataset, Dataset, ndarray[Any, dtype[int64]], NamedTuple]

Returns:

  • gains: GainTable Dataset containing the calibration solutions

  • model vis: Visibility Dataset containing model visibilities used in calibration

  • calibrated vis: Visibility Dataset containing calibrated visibilities

  • updated antenna flags: list of antenna flags used in final calibration step

  • calibration info: namedtuple of return variables:
    • corrcoeff: List of Pearson product-moment correlation coefficients for each polarisation of the calibrated visibilities and the model.

    • residual_std: List of residual visibility standard deviation per polarisation (K).

    • residual_max: List of residual visibility maximum absolution deviation per polarisation (K).

    • xy_phase: Estimated xy-phase error (degrees). Corrected for if

    • n_masked_initial: Initial number of masked antennas.

    • n_masked_final: Final number of masked antennas. jones_solve=False.

    • lst: Apparent sidereal time at station (degrees).

    • galactic_centre_elevation: Elevation of the galactic centre (degrees).

    • sun_elevation: Elevation of the sun (degrees).

    • sun_adjustment_factor: Solar adjustment factor (if adjust_solar_model=True, “1” otherwise).

ska_low_mccs_calibration.calibration.xy_phase_estimate(vis, modelvis)

Estimate uncalibrated xy-phase by comparing calibrated and model visibilities.

Parameters:
  • vis (Dataset) – Calibrated Visibility xarray

  • modelvis (Dataset) – Model Visibility xarray to update

Return type:

float

Returns:

xy-phase estimate (radians)

Embedded Element Pattern Support.

ska_low_mccs_calibration.eep.azel_to_eep_indices(eep_shape, eep_rotation_deg, azimuth_deg, elevation_deg)

Return nearest EEP pixels for a set of azimuths and elevations.

All elevations should be greater than or equal to 0. If the station is rotated in azimuth, the EEP pixel coordaintes are rotated accordingly.

Parameters:
  • eep_shape (ndarray, list) – shape of EEP azimuth and zenith angle axes, [Naz,Nza]

  • eep_rotation_deg (float) – rotation of the station (and therefore EEPs) in degrees

  • azimuth_deg (ndarray, list) – azimuth of pixels in degrees

  • elevation_deg (ndarray, list) – elevation of pixels in degrees

Return type:

list[ndarray[Any, dtype[int64]]]

Returns:

[azimuth indices, zenith angle indices]

ska_low_mccs_calibration.eep.convert_eep2npy(eep_filename, npy_dir='.', dtype=<class 'numpy.complex64'>)

Convert MATLAB EEP files to NumPy format and write them to disk.

EEPs are split into separate files for the phi and theta polarisations.

Parameters:
  • eep_filename (str) – MATLAB filename or list of names to load.

  • npy_dir (str) – Output directory for npy files. Default is the current working directory.

  • dtype (dtype) – The NumPy data-type to use for the output arrays. MATLAB files have dtype=complex128, but by default they will be reduced on output to dtype=complex64.

Return type:

None

ska_low_mccs_calibration.eep.get_nearest_file(filename_pattern, frequency_mhz)

Find the closest frequency of the available EEP files.

This method is very basic at this point and requires the filenames to contain the substring “<freq>MHz”.

Parameters:
  • filename_pattern (str) – Filename pattern with standard glob wildcards to search for closest frequency match. Filenames must contain sub-string “<freq>MHz”.

  • frequency_mhz (float) – Frequency of data being processed, in MHz.

Return type:

str

Returns:

Frequency of the nearest file.

ska_low_mccs_calibration.eep.load_eeps(frequency_mhz, path, filebase='FEKO_AAVS3_vogel_256_elem_50ohm_', dtype=<class 'numpy.complex64'>, conjugate=False)

Load EEP NumPy files into beam ndarrays of shape [az, za, nant].

Antennas are reordered to match the visibility data order, as specified by antenna_order. EEPs are stored in separate arrays for the four polarisation combinations: X,Y instrument polarisations and phi,theta sky polarisations.

Note: if EEP images are in MATLAB format, then before calling this function one should first call convert_eep2npy() to convert to NumPy format.

Parameters:
  • frequency_mhz (float) – Frequency of data being processed, in MHz.

  • path (str) – Directory path to EEP image files

  • filebase (str) – Start of filenames, before frequency specification. Files are assumed to have a specific form, which for polarisation combination X/phi is: path/filebase<freq>MHz_Xpol_phi.npy. Default is filebase=”FEKO_AAVS3_vogel_256_elem_50ohm_”.

  • antenna_order – The antenna order for these models in the visibility data. The first EEP image corresponds to antenna antenna_order[0] in the visibilties, for example. Default is antenna_order=None.

  • dtype (dtype) – The NumPy data-type to use to store the arrays. MATLAB files have dtype=complex128, but by default convert_eep2npy() reduces the npy output to dtype=complex64. Default here is dtype=complex64.

  • conjugate (bool) – Whether to conjugate the data. Default is False.

Return type:

list[ndarray[Any, dtype[Union[complex64, complex128]]]]

Returns:

[eep_x_phi, eep_x_theta, eep_y_phi, eep_y_theta]. One ndarray for each of the four polarisation combinations. Dimensions are [azimuth, zenith angle, antenna (in data order)].

ska_low_mccs_calibration.eep.reconfigure_eeps(eeps, eep_rotation_deg, azimuth_deg, elevation_deg)

Reconfigure EEPs for operations with the local sky model.

This amounts to a change of axis ordering and a normalisation with respect to pixel centres of the sky model. At present this is a simple nearest-pixel resampling.

On input the shape of the beam arrays is [azimuth, zenith angle, antenna]. On output the shape of the beam arrays is [antenna, zenith angle, azimuth]

Parameters:
  • eeps (list[ndarray[Any, dtype[Union[complex64, complex128]]]]) – list of four EEP sets: [Xphi, Xtheta, Yphi, Ytheta]

  • eep_rotation_deg (float) – rotation of the station (and therefore EEPs) in degrees

  • azimuth_deg (ndarray, list) – azimuth of new pixels in degrees

  • elevation_deg (ndarray, list) – elevation of new pixels in degrees

Return type:

list[ndarray[Any, dtype[Union[complex64, complex128]]]]

Returns:

updated beam arrays

ska_low_mccs_calibration.eep.reorder_eep_antennas(eeps, antenna_order)

Reorder antennas in EEPs of shape [az, za, nant].

Antennas are reordered to match the visibility data order, as specified by antenna_order. EEPs are stored in separate arrays for the four polarisation combinations: X,Y instrument polarisations and phi,theta sky polarisations.

Parameters:
  • eeps (list[ndarray[Any, dtype[Union[complex64, complex128]]]]) – list of four EEP sets: [Xphi, Xtheta, Yphi, Ytheta]

  • antenna_order (ndarray, list) – The antenna order for these models in the visibility data. The first EEP image corresponds to antenna antenna_order[0] in the visibilties, for example. Default is antenna_order=None.

Return type:

list[ndarray[Any, dtype[Union[complex64, complex128]]]]

Returns:

updated beam arrays

ska_low_mccs_calibration.eep.resample_eeps(eeps, eep_rotation_deg, azimuth_deg, elevation_deg)

Resample beam models at the LSM pixel centres.

At present this is a simple nearest-pixel resampling.

On input the shape of the beam arrays is [antenna, zenith angle, azimuth]. On output the shape of the beam arrays is [antenna, pixel]

Parameters:
  • eeps (list[ndarray[Any, dtype[Union[complex64, complex128]]]]) – list of four EEP sets: [Xphi, Xtheta, Yphi, Ytheta]

  • eep_rotation_deg (float) – rotation of the station (and therefore EEPs) in degrees

  • azimuth_deg (ndarray, list) – azimuth of new pixels in degrees

  • elevation_deg (ndarray, list) – elevation of new pixels in degrees

Return type:

list[ndarray[Any, dtype[Union[complex64, complex128]]]]

Returns:

updated beam arrays

Sky Model Support.

ska_low_mccs_calibration.sky_model.dft_sky_to_vis(eeps, index_array, ft_phasor, lsm)

Evaluate the DFT of the sky model for each baseline.

Parameters:
  • eeps (list[ndarray[Any, dtype[Union[complex64, complex128]]]]) – list of four EEP sets: [Xphi, Xtheta, Yphi, Ytheta]

  • index_array (ndarray, list) – first and second antenna indices for each baseline

  • ft_phasor (ndarray, list) – array of DFT phase terms between pixels and baselines

  • lsm (ndarray, list) – list of pixel values

Return type:

list[ndarray[Any, dtype[Union[complex64, complex128]]]]

Returns:

complex model visibilities [XX, XY, YX, YY]

ska_low_mccs_calibration.sky_model.generate_ft_phasor(delay_ns, frequency_mhz, azimuth_deg, elevation_deg)

Generate Fourier phase shifts for the LSM.

Genearte a complex array of phase shifts for the DFT from sky model pixels to sky model visibilities.

Parameters:
  • delay_ns (ndarray, list) – array of shape (3,nbaselines) containing the relative baseline delays in nanoseconds.

  • frequency_mhz (float) – observing frequency in MHz

  • azimuth_deg (ndarray, list) – array of shape (nskypixels) containing pixel azimuths.

  • elevation_deg (ndarray, list) – array of shape (nskypixels) containing zenith angles.

Return type:

ndarray[Any, dtype[Union[complex64, complex128]]]

Returns:

complex array of DFT phase terms, shape (nbaselines,nskypixels).

ska_low_mccs_calibration.sky_model.gsmap_lsm(gsmap, location, time, elevation_min=0.0)

Extract relevant pixels from the GSM.

Return the azimuth, elevation and temperature values of all pixels in the Healpix map above elevation_min degrees.

Parameters:
  • gsmap (ndarray, list) – Healpix map.

  • location (Quantity) – Astropy coordinates of the observation location.

  • time (Quantity) – Astropy time of the observation.

  • elevation_min (float) – min elevation cutoff in degrees (default is 0 deg)

Return type:

list[ndarray[Any, dtype[Union[float32, float64]]]]

Returns:

[azimuth_deg, elevation_deg, temperature]

ska_low_mccs_calibration.sky_model.gsmodel(frequency_mhz, nside)

Call PyGDSM to evaluate the diffuse galactic sky model.

PyGDSM model GlobalSkyModel16 is used.

Parameters:
Return type:

ndarray[Any, dtype[Union[float32, float64]]]

Returns:

List containing the Healpix map at each frequency

ska_low_mccs_calibration.sky_model.predict_vis(eeps, eep_rotation_deg, index_array, delay_ns, frequency_mhz, azimuth_deg, elevation_deg, lsm)

Evaluate the DFT of the sky model for each baseline.

Parameters:
  • eeps (list[ndarray[Any, dtype[Union[complex64, complex128]]]]) – list of four EEP sets: [Xphi, Xtheta, Yphi, Ytheta]

  • eep_rotation_deg (float) – rotation of the station (and therefore EEPs) in degrees

  • index_array (ndarray, list) – first and second antenna indices for each baseline

  • delay_ns (ndarray, list) – array of shape (3,nbaselines) containing the relative baseline delays in nanoseconds

  • frequency_mhz (float) – observing frequency in MHz

  • azimuth_deg (ndarray, list) – array of shape (nskypixels) containing pixel azimuths.

  • elevation_deg (ndarray, list) – array of shape (nskypixels) containing zenith angles.

  • lsm (ndarray, list) – list of pixel values

Return type:

list[ndarray[Any, dtype[Union[complex64, complex128]]]]

Returns:

complex model visibilities [XX, XY, YX, YY]

ska_low_mccs_calibration.sky_model.solar_lsm(frequency_mhz, location, time)

Get the position of the sun and estimate the temperature.

Temperature is based on the analytical expression for the quiet Sun in the range 30 - 350 MHz from Benz A. O. (2009) “Radio emission of the quiet sun.” Landolt Börnstein, 4B (103). This will be an underestimate when the sun is active.

Parameters:
  • frequency_mhz (ndarray, list, float) – Frequency or list of frequencies to evaluate (MHz).

  • location (Quantity) – Astropy coordinates of the observation location.

  • time (Quantity) – Astropy time of the observation.

Return type:

list[ndarray[Any, dtype[Union[float32, float64]]]]

Returns:

[azimuth_deg, elevation_deg, temperature]

Extra Utility Functions.

ska_low_mccs_calibration.utils.apply_gaintable(vis, gt, inverse=False)

Apply a GainTable to a Visibility.

This is a temporary local version of ska-sdp-func-python.operations.apply_gaintable that avoid a bug in the original. Will remove once the ska-sdp-func-python version is fixed. The bug is that the function does not transpose the rightmost matrix. This may not always be the desired action, for instance if calibration was done separately for each polarised visibility type, but is required when the solutions are Jones matrices.

Note: this is a temporary function and has not been made to robustly handle all situations. For instance, it ignores flags and does not properly handle sub-bands and partial solution intervals. It also assumes that the matrices being applied are Jones matrices and the general ska-sdp-datamodels xarrays are used.

Parameters:
  • vis (Dataset) – ska-sdp-datamodels Visibility xarray containing visibilities.

  • gt (Dataset) – ska-sdp-datamodels GainTable xarray containing gain solutions.

  • inverse (bool) – Apply the inverse of the gain. This requires the gain matrices to be square. (default=False)

Return type:

Dataset

Returns:

ska-sdp-datamodels Visibility xarray containing calibrated visibilities.

ska_low_mccs_calibration.utils.enu_to_geocentric(antpos_enu, location)

Convert East,North,Up antenna coordinates to geocentric.

Convert antenna coordinates from the local frame to geocentric.

Parameters:
  • antpos_enu (ndarray, list) – array of East,North,Up coordinates. This array must have shape […, 3], with East,North,Up coordinates along the last axis.

  • location (Quantity) – Astropy coordinates of the observation location.

Return type:

ndarray[Any, dtype[Union[float32, float64]]]

Returns:

ndarray with the same shape as antpos_enu and X,Y,Z coordinates along the last axis.

ska_low_mccs_calibration.utils.read_yaml_config(yaml_file)

Extract configuration metadata from observation YAML file.

Parameters:

yaml_file (str) – YAML file to read

Return type:

tuple[Quantity, ndarray[Any, dtype[bool_]], tuple[ndarray[Any, dtype[int64]], ndarray[Any, dtype[int64]]], tuple[ndarray[Any, dtype[Union[float32, float64]]], ndarray[Any, dtype[Union[float32, float64]]], ndarray[Any, dtype[Union[float32, float64]]]], list[int], float]

Returns:

tuple of configuration parameters: location, # Astropy array location antenna_masks, # bool list of flagged antennas [vis ant1, vis ant 2], # first and second antenna for each baseline [E, N, U], # Local East, North, Up location of each antenna [XX, XY, YX, YY], # index of each polarisation in vis dataset rotation, # azimuthal rotation of station (bearing clockwise from north)

ska_low_mccs_calibration.utils.sdp_visibility_datamodel(vis, flags, uvw, ant1, ant2, location, antpos_enu, time, int_time, frequency_mhz, channel_bandwidth_mhz=0.78125)

Package MCCS data in a ska-sdp-datamodels Visibility xarray.

Currently set up for a single time step and frequency channel

Parameters:
  • vis (ndarray[Any, dtype[Union[complex64, complex128]]]) – 2D or 3D array of visibilties, shape [<ntimes>,nbaselines,npol

  • flags (ndarray[Any, dtype[bool_]]) – 2D or 3D array of bad-data flags, shape [<ntimes>,nbaselines,npol]

  • uvw (ndarray[Any, dtype[Union[float32, float64]]]) – 2D or 3D array of antenna coordinates in metres, shape [<ntimes>,nbaselines,3]

  • ant1 (ndarray, list) – 1D array containing the indices of the first antenna in each baseline

  • ant2 (ndarray, list) – 1D array containing the indices of the second antenna in each baseline

  • location (Quantity) – Astropy coordinates of the observation location.

  • antpos_enu (ndarray, list) – Array of antenna East,North,Up coordinates in metres. This array must have shape [nant, 3], with East,North,Up coordinates along the 2nd axis.

  • time (Quantity) – Astropy time of the observation time step.

  • int_time (float) – Integration time of data samples (seconds).

  • frequency_mhz (float) – Observing frequency of data channel (MHz).

  • channel_bw_MHz – Channel bandwidth (MHz). Default = 400MHz/512 = 0.78125

Return type:

Dataset

Returns:

ska-sdp-datamodels Visibility xarray containing visibilities and metadata