fixms.fix_ms_corrs#
Fix the correlation rotation of ASKAP MSs.
Converts the ASKAP standard correlations to the ‘standard’ correlations This will make them compatible with most imagers (e.g. wsclean, CASA)
The new correlations are placed in a new column called ‘CORRECTED_DATA’
Attributes#
Functions#
|
Check if the data in the data_column, with a correction applied, matches |
|
Check if the data in the data_column, with a correction applied, matches |
|
|
|
Convert ASKAP standard correlations to the 'standard' correlations |
|
Apply corrections to the ASKAP visibilities to bring them inline with |
|
Apply corrections to the ASKAP visibilities to bring them inline with |
|
Gather with a limit on the number of coroutines running at once. |
Generator function that will yield a chunk of data from the ms data table. |
|
|
Returns the number of chunks that are needed to iterator over the datacolumsn |
|
Get the polarization axis from the ASKAP MS. Checks are performed |
|
Get the polarization axis from the ASKAP MS. Checks are performed |
|
Work on a chunk of data and write it back to the MS |
|
Module Contents#
- fixms.fix_ms_corrs.check_data(ms: pathlib.Path, data_column: str, corrected_data_column: str) bool[source]#
Check if the data in the data_column, with a correction applied, matches the data in the corrected_data_column.
- Parameters:
ms (Path) – MeasurementSet to check
data_column (str) – Data column to check
corrected_data_column (str) – Corrected data column to check
- Returns:
If the data matches
- Return type:
bool
- async fixms.fix_ms_corrs.check_data_coro(ms: pathlib.Path, data_column: str, corrected_data_column: str) Awaitable[bool][source]#
Check if the data in the data_column, with a correction applied, matches the data in the corrected_data_column.
- Parameters:
ms (Path) – MeasurementSet to check
data_column (str) – Data column to check
corrected_data_column (str) – Corrected data column to check
- Returns:
If the data matches
- Return type:
bool
- fixms.fix_ms_corrs.convert_correlations(correlations: numpy.ndarray, pol_axis: astropy.units.Quantity, fix_stokes_factor: bool = True) numpy.ndarray[source]#
Convert ASKAP standard correlations to the ‘standard’ correlations
- Parameters:
correlations (np.ndarray) – The correlations from the MS. Has shape (NCOR, NCHAN, 4)
pol_axis (astropy.units.Quantity) – The polarization axis angle of the MS
fix_stokes_factor (bool, optional) – Whether to fix the Stokes factor. Defaults to True.
- Returns:
The correlations in the ‘standard’ format
- Return type:
np.ndarray
NOTES: In general, ASKAP forms Stokes I, Q, U, V as follows:
⎡I⎤ ⎡ 1 0 0 1 ⎤ ⎡XXₐ⎤ ⎢ ⎥ ⎢ ⎥ ⎢ ⎥ ⎢Q⎥ ⎢sin(2⋅P) cos(2⋅P) cos(2⋅P) -sin(2⋅P)⎥ ⎢XYₐ⎥ ⎢ ⎥ = ⎢ ⎥⋅⎢ ⎥ ⎢U⎥ ⎢-cos(2⋅P) sin(2⋅P) sin(2⋅P) cos(2⋅P) ⎥ ⎢YXₐ⎥ ⎢ ⎥ ⎢ ⎥ ⎢ ⎥ ⎣V⎦ ⎣ 0 -1.0⋅ⅈ 1.0⋅ⅈ 0 ⎦ ⎣YYₐ⎦
Where P is the polarization axis angle. In the common case of P=-45deg this becomes:
⎡I⎤ ⎡1 0 0 1⎤ ⎡XXₐ⎤ ⎢ ⎥ ⎢ ⎥ ⎢ ⎥ ⎢Q⎥ ⎢-1 0 0 1⎥ ⎢XYₐ⎥ ⎢ ⎥ = ⎢ ⎥⋅⎢ ⎥ ⎢U⎥ ⎢0 -1 -1 0⎥ ⎢YXₐ⎥ ⎢ ⎥ ⎢ ⎥ ⎢ ⎥ ⎣V⎦ ⎣0 -1.0⋅ⅈ 1.0⋅ⅈ 0⎦ ⎣YYₐ⎦
or
⎡I⎤ ⎡ XXₐ + YYₐ ⎤ ⎢ ⎥ ⎢ ⎥ ⎢Q⎥ ⎢ -XXₐ + YYₐ ⎥ ⎢ ⎥ = ⎢ ⎥ ⎢U⎥ ⎢ -XYₐ - YXₐ ⎥ ⎢ ⎥ ⎢ ⎥ ⎣V⎦ ⎣-ⅈ⋅XYₐ + 1.0⋅ⅈ⋅YXₐ⎦
If instead the P=+45deg then this becomes:
⎡I⎤ ⎡1 0 0 1 ⎤ ⎡XXₐ⎤ ⎢ ⎥ ⎢ ⎥ ⎢ ⎥ ⎢Q⎥ ⎢1 0 0 -1⎥ ⎢XYₐ⎥ ⎢ ⎥ = ⎢ ⎥⋅⎢ ⎥ ⎢U⎥ ⎢0 1 1 0 ⎥ ⎢YXₐ⎥ ⎢ ⎥ ⎢ ⎥ ⎢ ⎥ ⎣V⎦ ⎣0 -1.0⋅ⅈ 1.0⋅ⅈ 0 ⎦ ⎣YYₐ⎦
or
⎡I⎤ ⎡ XXₐ + YYₐ ⎤ ⎢ ⎥ ⎢ ⎥ ⎢Q⎥ ⎢ XXₐ - YYₐ ⎥ ⎢ ⎥ = ⎢ ⎥ ⎢U⎥ ⎢ XYₐ + YXₐ ⎥ ⎢ ⎥ ⎢ ⎥ ⎣V⎦ ⎣-ⅈ⋅XYₐ + 1.0⋅ⅈ⋅YXₐ⎦
However, most imagers (e.g. wsclean, CASA) expect
⎡I⎤ ⎡0.5 0 0 0.5 ⎤ ⎡XX_w⎤ ⎢ ⎥ ⎢ ⎥ ⎢ ⎥ ⎢Q⎥ ⎢0.5 0 0 -0.5⎥ ⎢XY_w⎥ ⎢ ⎥ = ⎢ ⎥⋅⎢ ⎥ ⎢U⎥ ⎢ 0 0.5 0.5 0 ⎥ ⎢YX_w⎥ ⎢ ⎥ ⎢ ⎥ ⎢ ⎥ ⎣V⎦ ⎣ 0 -0.5⋅i 0.5⋅i 0 ⎦ ⎣YY_w⎦
or
⎡I⎤ ⎡ 0.5⋅XX_w + 0.5⋅YY_w ⎤ ⎢ ⎥ ⎢ ⎥ ⎢Q⎥ ⎢ 0.5⋅XX_w - 0.5⋅YY_w ⎥ ⎢ ⎥ = ⎢ ⎥ ⎢U⎥ ⎢ 0.5⋅XY_w + 0.5⋅YX_w ⎥ ⎢ ⎥ ⎢ ⎥ ⎣V⎦ ⎣-0.5⋅i⋅XY_w + 0.5⋅i⋅YX_w⎦
To convert between the two, we can use the following matrix:
⎡XX_w⎤ ⎡sin(2.0⋅P) + 1 cos(2.0⋅P) cos(2.0⋅P) 1 - sin(2.0⋅P)⎤ ⎡XXₐ⎤ ⎢ ⎥ ⎢ ⎥ ⎢ ⎥ ⎢XY_w⎥ ⎢ -cos(2.0⋅P) sin(2.0⋅P) + 1 sin(2.0⋅P) - 1 cos(2.0⋅P) ⎥ ⎢XYₐ⎥ ⎢ ⎥ = ⎢ ⎥⋅⎢ ⎥ ⎢YX_w⎥ ⎢ -cos(2.0⋅P) sin(2.0⋅P) - 1 sin(2.0⋅P) + 1 cos(2.0⋅P) ⎥ ⎢YXₐ⎥ ⎢ ⎥ ⎢ ⎥ ⎢ ⎥ ⎣YY_w⎦ ⎣1 - sin(2.0⋅P) -cos(2.0⋅P) -cos(2.0⋅P) sin(2.0⋅P) + 1⎦ ⎣YYₐ⎦
Where _w is the ‘wsclean’ format and _a is the ‘ASKAP’ format.
In the case of P=-45deg this becomes:
⎡XX_w⎤ ⎡0 0 0 2⎤ ⎡XXₐ⎤ ⎢ ⎥ ⎢ ⎥ ⎢ ⎥ ⎢XY_w⎥ ⎢0 0 -2 0⎥ ⎢XYₐ⎥ ⎢ ⎥ = ⎢ ⎥⋅⎢ ⎥ ⎢YX_w⎥ ⎢0 -2 0 0⎥ ⎢YXₐ⎥ ⎢ ⎥ ⎢ ⎥ ⎢ ⎥ ⎣YY_w⎦ ⎣2 0 0 0⎦ ⎣YYₐ⎦
or
⎡XX_w⎤ ⎡2⋅YYₐ ⎤ ⎢ ⎥ ⎢ ⎥ ⎢XY_w⎥ ⎢-2⋅YXₐ⎥ ⎢ ⎥ = ⎢ ⎥ ⎢YX_w⎥ ⎢-2⋅XYₐ⎥ ⎢ ⎥ ⎢ ⎥ ⎣YY_w⎦ ⎣2⋅XXₐ ⎦
And in the case of P=+45deg this becomes:
⎡XX_w⎤ ⎡2 0 0 0⎤ ⎡XXₐ⎤ ⎢ ⎥ ⎢ ⎥ ⎢ ⎥ ⎢XY_w⎥ ⎢0 2 0 0⎥ ⎢XYₐ⎥ ⎢ ⎥ = ⎢ ⎥⋅⎢ ⎥ ⎢YX_w⎥ ⎢0 0 2 0⎥ ⎢YXₐ⎥ ⎢ ⎥ ⎢ ⎥ ⎢ ⎥ ⎣YY_w⎦ ⎣0 0 0 2⎦ ⎣YYₐ⎦
or
⎡XX_w⎤ ⎡2⋅XXₐ⎤ ⎢ ⎥ ⎢ ⎥ ⎢XY_w⎥ ⎢2⋅XYₐ⎥ ⎢ ⎥ = ⎢ ⎥ ⎢YX_w⎥ ⎢2⋅YXₐ⎥ ⎢ ⎥ ⎢ ⎥ ⎣YY_w⎦ ⎣2⋅YYₐ⎦
That is, when P=+45 H and V are perfectly aligned with X and Y, and only the factor of 2 is required.
- fixms.fix_ms_corrs.fix_ms_corrs(ms: pathlib.Path, chunksize: int = 1000, max_chunks: int = 1000, data_column: str = 'DATA', corrected_data_column: str = 'CORRECTED_DATA', fix_stokes_factor: bool = True) None[source]#
Apply corrections to the ASKAP visibilities to bring them inline with what is expectede from other imagers, including CASA and WSClean. The original data in data_column are copied to corrected_data_column before the correction is applied. This is done to ensure that the original data are not lost.
If ‘corrected_data_column’ is detected as an existing column then the correction will not be applied.
- Parameters:
ms (Path) – Path of the ASKAP measurement set tto correct.
chunksize (int, optional) – Size of chunked data to correct. Defaults to 10_000.
data_column (str, optional) – The name of the data column to correct. Defaults to “DATA”.
corrected_data_column (str, optional) – The name of the corrected data column. Defaults to “CORRECTED_DATA”.
fix_stokes_factor (bool, optional) – Whether to fix the Stokes factor. Defaults to True.
- async fixms.fix_ms_corrs.fix_ms_corrs_coro(ms: pathlib.Path, chunksize: int = 1000, max_chunks: int = 1000, data_column: str = 'DATA', corrected_data_column: str = 'CORRECTED_DATA', fix_stokes_factor: bool = True) Awaitable[source]#
Apply corrections to the ASKAP visibilities to bring them inline with what is expectede from other imagers, including CASA and WSClean. The original data in data_column are copied to corrected_data_column before the correction is applied. This is done to ensure that the original data are not lost.
If ‘corrected_data_column’ is detected as an existing column then the correction will not be applied.
- Parameters:
ms (Path) – Path of the ASKAP measurement set tto correct.
chunksize (int, optional) – Size of chunked data to correct. Defaults to 10_000.
data_column (str, optional) – The name of the data column to correct. Defaults to “DATA”.
corrected_data_column (str, optional) – The name of the corrected data column. Defaults to “CORRECTED_DATA”.
fix_stokes_factor (bool, optional) – Whether to fix the Stokes factor. Defaults to True.
- async fixms.fix_ms_corrs.gather_with_limit(limit: int, *coros: Awaitable)[source]#
Gather with a limit on the number of coroutines running at once.
- Parameters:
limit (int) – The number of coroutines to run at once
coros (Awaitable) – The coroutines to run
- Returns:
The result of the coroutines
- Return type:
Awaitable
- async fixms.fix_ms_corrs.get_data_chunk_generator(ms: pathlib.Path, chunksize: int, data_column: str = 'DATA') AsyncGenerator[numpy.ndarray, None][source]#
Generator function that will yield a chunk of data from the ms data table.
- Parameters:
ms (Path) – Measurement set whose data will be iterated over
chunksize (int) – The number of rows to process per chunk
data_column (str, optional) – The column name of the data to iterate. Defaults to “DATA”.
- Yields:
Iterator[np.ndarray] – Chunk of datta to process
- async fixms.fix_ms_corrs.get_nchunks_coro(ms: pathlib.Path, chunksize: int, data_column: str = 'DATA') Awaitable[int][source]#
Returns the number of chunks that are needed to iterator over the datacolumsn using a specified chunksize.
- Parameters:
ms (Path) – Measurement sett thatt will be iterated over
chunksize (int) – Size of a single chunk
data_column (str, optional) – Name of the datacolumn that will be iterated over. Defaults to “DATA”.
- Returns:
Number of chunks in a measurement set
- Return type:
int
- fixms.fix_ms_corrs.get_pol_axis(ms: pathlib.Path, feed_idx: int | None = None, col='RECEPTOR_ANGLE') astropy.units.Quantity[source]#
Get the polarization axis from the ASKAP MS. Checks are performed to ensure this polarisation axis angle is constant throughout the observation.
- Parameters:
ms (Path) – The path to the measurement set that will be inspected
feed_idx (Optional[int], optional) – Specify the entry in the FEED
table of `ms` to return. This might be required when a subset of a
measurement set has been extracted from an observation with a varying
orientation.
- Returns:
The rotation of the PAF throughout the observing.
- Return type:
astropy.units.Quantity
- async fixms.fix_ms_corrs.get_pol_axis_coro(ms: pathlib.Path, feed_idx: int | None = None, col='RECEPTOR_ANGLE') Awaitable[astropy.units.Quantity][source]#
Get the polarization axis from the ASKAP MS. Checks are performed to ensure this polarisation axis angle is constant throughout the observation.
- Parameters:
ms (Path) – The path to the measurement set that will be inspected
feed_idx (Optional[int], optional) – Specify the entry in the FEED
table of `ms` to return. This might be required when a subset of a
measurement set has been extracted from an observation with a varying
orientation.
- Returns:
The rotation of the PAF throughout the observing.
- Return type:
astropy.units.Quantity
- async fixms.fix_ms_corrs.process_chunk(tab: casacore.tables.table, data_column: str, corrected_data_column: str, pol_axis: astropy.units.Quantity, fix_stokes_factor: bool, start_row: int, chunksize: int, pbar: tqdm.auto.tqdm, chunk_num: int)[source]#
Work on a chunk of data and write it back to the MS
- Parameters:
data_chunk (np.ndarray) – Chunk of MS rows
pol_axis (u.Quantity) – Rotation angle of the polarization axis
fix_stokes_factor (bool) – Whether to fix the Stokes factor
tab (table) – MS table object
corrected_data_column (str) – Corrected data column
start_row (int) – Starting row of the chunk