Data¶
Data¶
Memory Data models¶
The data models used in ARL:
Note
There are two visibility formats:
BlockVisibility
is conceived as an ingest and calibration format. The visibility
data are kept in a block of shape (number antennas, number antennas, number channels,
number polarisation). One block is kept per integration. The other columns are time and uvw.
The sampling in time is therefore the same for all baselines.
Visibility
is designed to hold coalesced data where the integration time and
channel width can vary with baseline length. The visibility data are kept in a visibility
vector of length equal to the number of polarisations. Everything else is a separate
column: time, frequency, uvw, channel_bandwidth, integration time.
-
class
BlockVisibility
(data=None, frequency=None, channel_bandwidth=None, phasecentre=None, configuration=None, uvw=None, time=None, vis=None, weight=None, integration_time=None, polarisation_frame=<data_models.polarisation.PolarisationFrame object>)[source]¶ Block Visibility table class
BlockVisibility with uvw, time, integration_time, frequency, channel_bandwidth, pol, a1, a2, vis, weight Columns in a numpy structured array.
BlockVisibility is defined to hold an observation with one direction.
The phasecentre is the direct of delay tracking i.e. n=0. If uvw are rotated then this should be updated with the new delay tracking centre. This is important for wstack and wproject algorithms.
Polarisation frame is the same for the entire data set and can be stokesI, circular, linear
The configuration is also an attribute
-
class
Configuration
(name='', data=None, location=None, names='%s', xyz=None, mount='alt-az', frame=None, receptor_frame=<data_models.polarisation.ReceptorFrame object>, diameter=None)[source]¶ Describe a Configuration as locations in x,y,z, mount type, diameter, names, and overall location
-
diameter
¶ diameter of antennas/stations
-
mount
¶ Mount type
-
names
¶ Names of the antennas/stations
-
xyz
¶ XYZ locations of antennas/stations
-
-
class
GainTable
(data=None, gain: <built-in function array> = None, time: <built-in function array> = None, interval=None, weight: <built-in function array> = None, residual: <built-in function array> = None, frequency: <built-in function array> = None, receptor_frame: data_models.polarisation.ReceptorFrame = <data_models.polarisation.ReceptorFrame object>)[source]¶ Gain table with data_models: time, antenna, gain[:, chan, rec, rec], weight columns
The weight is usually that output from gain solvers.
-
class
Image
[source]¶ Image class with Image data (as a numpy.array) and the AstroPy implementation of a World Coodinate System
Many operations can be done conveniently using numpy libs on Image.data_models.
Most of the imaging libs require an image in canonical format: - 4 axes: RA, DEC, POL, FREQ
The conventions for indexing in WCS and numpy are opposite. - In astropy.wcs, the order is (longitude, latitude, polarisation, frequency) - in numpy, the order is (frequency, polarisation, latitude, longitude)
Warning
The polarisation_frame is kept in two places, the WCS and the polarisation_frame variable. The latter should be considered definitive.
-
class
Skycomponent
(direction=None, frequency=None, name=None, flux=None, shape='Point', polarisation_frame=<data_models.polarisation.PolarisationFrame object>, params=None)[source]¶ Skycomponents are used to represent compact sources on the sky. They possess direction, flux as a function of frequency and polarisation, shape (with params), and polarisation frame.
For example, the following creates and predicts the visibility from a collection of point sources drawn from the GLEAM catalog:
sc = create_low_test_skycomponents_from_gleam(flux_limit=1.0, polarisation_frame=PolarisationFrame("stokesIQUV"), frequency=frequency, kind='cubic', phasecentre=phasecentre, radius=0.1) model = create_image_from_visibility(vis, cellsize=0.001, npixel=512, frequency=frequency, polarisation_frame=PolarisationFrame('stokesIQUV')) bm = create_low_test_beam(model=model) sc = apply_beam_to_skycomponent(sc, bm) vis = predict_skycomponent_visibility(vis, sc)
-
class
Visibility
(data=None, frequency=None, channel_bandwidth=None, phasecentre=None, configuration=None, uvw=None, time=None, antenna1=None, antenna2=None, vis=None, weight=None, imaging_weight=None, integration_time=None, polarisation_frame=<data_models.polarisation.PolarisationFrame object>, cindex=None, blockvis=None)[source]¶ Visibility table class
Visibility with uvw, time, integration_time, frequency, channel_bandwidth, a1, a2, vis, weight as separate columns in a numpy structured array, The fundemental unit is a complex vector of polarisation.
Visibility is defined to hold an observation with one direction. Polarisation frame is the same for the entire data set and can be stokes, circular, linear The configuration is also an attribute
The phasecentre is the direct of delay tracking i.e. n=0. If uvw are rotated then this should be updated with the new delay tracking centre. This is important for wstack and wproject algorithms.
If a visibility is created by coalescence then the cindex column is filled with a pointer to the row in the original block visibility that this row has a value for. The original blockvisibility is also preserves as n attribute so that decoalescence is expedited. If you don’t need that then the storage can be released by setting self.blockvis to None
Buffer Data Models¶
Buffer equivalents of memory data models
To create from an existing model in the buffer, make some changes, and then sync to buffer:
my_buffer_skymodel = BufferSkyModel(conf["buffer"], conf["inputs"]["skymodel"])
my_memory_skymodel = my_buffer_skymodel.memory_data_model
... do some stuff
my_memory_skymodel.sync()
To create a new buffer for a memory_data_model:
my_buffer_skymodel = BufferSkyModel(conf["buffer"], conf["inputs"]["skymodel"], my_memory_skymodel)
my_memory_skymodel.sync()
An explicit sync is required in both cases
-
class
BufferBlockVisibility
(json_buffer, json_model, mdm=None)[source]¶ Buffer version of memory data model BlockVisibility
-
class
BufferDataModel
(json_buffer, json_model, mdm=None)[source]¶ Buffer version of data model
To create from an existing model in the buffer, make some changes, and then sync to buffer:
my_buffer_skymodel = BufferSkyModel(conf["buffer"], conf["inputs"]["skymodel"]) my_memory_skymodel = my_buffer_skymodel.memory_data_model() ... do some stuff my_memory_skymodel.sync()
To create a new buffer for a memory_data_model:
my_buffer_skymodel = BufferSkyModel(conf["buffer"], conf["inputs"]["skymodel"], my_memory_skymodel) my_memory_skymodel.sync()
An explicit sync is required in both cases.
-
class
BufferGainTable
(json_buffer, json_model, mdm=None)[source]¶ Buffer version of memory data model GainTable
Data model persistence¶
Functions to help with persistence of data models
These do data conversion and persistence. Functions from libs and processing_components are used.
-
buffer_data_model_to_memory
(jbuff, dm)[source]¶ Copy a buffer data model into memory data model
The file type is derived from the file extension. All are hdf only with the exception of Imaghe which can also be fits.
Parameters: - jbuff – JSON describing buffer
- dm – JSON describing data model
Returns: data model
-
convert_blockvisibility_to_hdf
(vis: data_models.memory_data_models.BlockVisibility, f)[source]¶ Convert blockvisibility to HDF
Parameters: - vis –
- f – HDF root
Returns:
-
convert_configuration_to_hdf
(config: data_models.memory_data_models.Configuration, f)[source]¶ Parameters: - config –
- f –
Returns:
-
convert_direction_from_string
(s: str)[source]¶ Convert direction (SkyCoord) from string
TODO: Make more general!
Parameters: s – String Returns:
-
convert_direction_to_string
(d: astropy.coordinates.sky_coordinate.SkyCoord)[source]¶ Convert SkyCoord to string
TODO: Make more general!
Parameters: d – SkyCoord Returns:
-
convert_earthlocation_from_string
(s: str)[source]¶ Convert Earth Location to string
Parameters: s – String Returns:
-
convert_earthlocation_to_string
(el: astropy.coordinates.earth.EarthLocation)[source]¶ Convert Earth Location to string
Parameters: el – Returns:
-
convert_gaintable_to_hdf
(gt: data_models.memory_data_models.GainTable, f)[source]¶ Convert GainTable to HDF
Parameters: - gt –
- f – HDF root
Returns:
-
convert_hdf_to_blockvisibility
(f)[source]¶ Convert HDF root to blockvisibility
Parameters: f – Returns:
-
convert_image_to_hdf
(im: data_models.memory_data_models.Image, f)[source]¶ Convert Image to HDF
Parameters: - im – Image
- f – HDF root
Returns:
-
convert_skycomponent_to_hdf
(sc: data_models.memory_data_models.Skycomponent, f)[source]¶ Convert Skycomponent to HDF :param sc: SkyComponent :param f: HDF root :return:
-
convert_visibility_to_hdf
(vis, f)[source]¶ Convert visibility to HDF
Parameters: - vis –
- f – HDF root
Returns:
-
export_blockvisibility_to_hdf5
(vis, filename)[source]¶ Export a BlockVisibility to HDF5 format
Parameters: - vis –
- filename –
Returns:
-
export_gaintable_to_hdf5
(gt: data_models.memory_data_models.GainTable, filename)[source]¶ Export a GainTable to HDF5 format
Parameters: - gt –
- filename –
Returns:
-
export_image_to_hdf5
(im, filename)[source]¶ Export an Image to HDF5 format
Parameters: - im –
- filename –
Returns:
-
export_skycomponent_to_hdf5
(sc: data_models.memory_data_models.Skycomponent, filename)[source]¶ Export a Skycomponent to HDF5 format
Parameters: - sc – SkyComponent
- filename –
Returns:
-
export_skymodel_to_hdf5
(sm, filename)[source]¶ Export a Skymodel to HDF5 format
Parameters: - sm –
- filename –
Returns:
-
export_visibility_to_hdf5
(vis, filename)[source]¶ Export a Visibility to HDF5 format
Parameters: - vis –
- filename –
Returns:
-
import_blockvisibility_from_hdf5
(filename)[source]¶ Import a Visibility from HDF5 format
Parameters: filename – Returns: If only one then a BlockVisibility, otherwise a list of BlockVisibility’s
-
import_gaintable_from_hdf5
(filename)[source]¶ Import GainTable(s) from HDF5 format
Parameters: filename – Returns: single gaintable or list of gaintables
-
import_image_from_hdf5
(filename)[source]¶ Import Image(s) from HDF5 format
Parameters: filename – Returns: single image or list of images
-
import_skycomponent_from_hdf5
(filename)[source]¶ Import Skycomponent(s) from HDF5 format
Parameters: filename – Returns: single skycomponent or list of skycomponents
-
import_skymodel_from_hdf5
(filename)[source]¶ Import a Skymodel from HDF5 format
Parameters: filename – Returns: SkyModel
-
import_visibility_from_hdf5
(filename)[source]¶ Import a Visibility from HDF5 format
Parameters: filename – Returns: If only one then a Visibility, otherwise a list of Visibilitys
-
memory_data_model_to_buffer
(model, jbuff, dm)[source]¶ Copy a memory data model to a buffer data model
The file type is derived from the file extension. All are hdf only with the exception of Imaghe which can also be fits.
Parameters: - model – Memory data model to be sent to buffer
- jbuff – JSON describing buffer
- dm – JSON describing data model
Parameter handling¶
We use the standard kwargs mechanism for arguments. For example:
kernelname = get_parameter(kwargs, "kernel", "2d")
oversampling = get_parameter(kwargs, "oversampling", 8)
padding = get_parameter(kwargs, "padding", 2)
The kwargs may need to be passed down to called functions.
All functions possess an API which is always of the form:
def processing_function(idatastruct1, idatastruct2, ..., *kwargs):
return odatastruct1, odatastruct2,... other
Processing parameters are passed via the standard Python kwargs approach.
Inside a function, the values are retrieved can be accessed directly from the kwargs dictionary, or if a default is needed a function can be used:
log = get_parameter(kwargs, 'log', None)
vis = get_parameter(kwargs, 'visibility', None)
Function parameters should obey a consistent naming convention:
Name | Meaning |
---|---|
vis | Name of Visibility |
sc | Name of Skycomponent |
gt | Name of GainTable |
conf | Name of Configuration |
im | Name of input image |
qa | Name of quality assessment |
log | Name of processing log |
If a function argument has a better, more descriptive name e.g. normalised_gt, newphasecentre, use it.
Keyword=value pairs should have descriptive names. The names should be lower case with underscores to separate words:
Name | Meaning | Example |
---|---|---|
loop_gain | Clean loop gain | 0.1 |
niter | Number of iterations | 10000 |
eps | Fractional tolerance | 1e-6 |
threshold | Absolute threshold | 0.001 |
fractional_threshold | Threshold as fraction of e.g. peak | 0.1 |
G_solution_interval | Solution interval for G term | 100 |
phaseonly | Do phase-only solutions | True |
phasecentre | Phase centre (usually as SkyCoord) | SkyCoord(“-1.0d”, “37.0d”, frame=’icrs’, equinox=’J2000’) |
spectral_mode | Visibility processing mode | ‘mfs’ or ‘channel’ |
Polarisation¶
- Functions for defining polarisation conventions. These include definitions via classes and
- conversion functions.
For example:
stokes = numpy.array(random.uniform(-1.0, 1.0, [3, 4, 128, 128]))
ipf = PolarisationFrame('stokesIQUV')
opf = PolarisationFrame('circular')
cir = convert_pol_frame(stokes, ipf, opf)
st = convert_pol_frame(cir, opf, ipf)
or:
stokes = numpy.array([1, 0.5, 0.2, -0.1])
circular = convert_stokes_to_circular(stokes)
These function operate on Numpy arrays. These are packaged for use in Images. The Image functions are probably more useful.
-
class
PolarisationFrame
(name)[source]¶ Define polarisation frames post correlation
-
npol
¶ Number of correlated polarisations
-
-
class
ReceptorFrame
(name)[source]¶ Define polarisation frames for receptors
circular, linear, and stokesI. The latter is non-physical but useful for some types of testing.
-
nrec
¶ Number of receptors (should be 2)
-
-
congruent_polarisation
(rec_frame: data_models.polarisation.ReceptorFrame, polarisation_frame: data_models.polarisation.PolarisationFrame)[source]¶ Are these receptor and polarisation frames congruent?
-
convert_circular_to_stokes
(circular, polaxis=1)[source]¶ Convert Circular to Stokes IQUV
Parameters: - circular – […,4] linear vector in RR, RL, LR, LL sequence
- polaxis – Axis of circular with polarisation (default 1)
Returns: Complex I,Q,U,V
Equation 4.58 TMS, inverted with numpy.linalg.inv
-
convert_linear_to_stokes
(linear, polaxis=1)[source]¶ Convert Linear to Stokes IQUV
Parameters: - linear – […,4] linear vector in XX, XY, YX, YY sequence
- polaxis – Axis of linear with polarisation (default 1)
Returns: Complex I,Q,U,V
Equation 4.58 TMS, inverted with numpy.linalg.inv
-
convert_stokes_to_circular
(stokes, polaxis=1)[source]¶ Convert Stokes IQUV to Circular
Parameters: - stokes – […,4] Stokes vector in I,Q,U,V (can be complex)
- polaxis – Axis of stokes with polarisation (default 1)
Returns: circular vector in RR, RL, LR, LL sequence
Equation 4.59 TMS
-
convert_stokes_to_linear
(stokes, polaxis=1)[source]¶ Convert Stokes IQUV to Linear
Parameters: - stokes – […,4] Stokes vector in I,Q,U,V (can be complex)
- polaxis – Axis of stokes with polarisation (default 1)
Returns: linear vector in XX, XY, YX, YY sequence
Equation 4.58 TMS
-
correlate_polarisation
(rec_frame: data_models.polarisation.ReceptorFrame)[source]¶ Gives the polarisation frame corresponding to a receptor frame
Parameters: rec_frame – Receptor frame Returns: PolarisationFrame
-
polmatrixmultiply
(cm, vec, polaxis=1)[source]¶ Matrix multiply of appropriate axis of vec […,:] by cm
For an image vec has axes [nchan, npol, ny, nx] and polaxis=1 For visibility vec has axes [row, nchan, npol] and polaxis=2
Parameters: - cm – matrix to apply
- vec – array to be multiplied […,:]
- polaxis – which axis contains the polarisation
Returns: multiplied vec