Calibration¶
This module provides basic radio interferometry calibration utilities. Calibration is the process of estimating the \(2\times 2\) Jones matrices which describe transformations of the signal as it propagates from source to observer. Currently, all utilities assume a discretised form of the radio interferometer measurement equation (RIME) as described in Radio Interferometer Measurement Equation.
Calibration is usually divided into three phases viz.
- First generation calibration (1GC): using an external calibrator to infer the gains during the target observation. Sometimes also refered to as calibrator transfer
- Second generation calibration (2GC): using a partially incomplete sky model to perform direction independent calibration. Also known as direction independent self-calibration.
- Third generation calibration (3GC): using a partially incomplete sky model to perform direction dependent calibration. Also known as direction dependent self-calibration.
On top of these three phases, there are usually
three possible calibration scenarios. The first
is when both the Jones terms and the visibilities
are assumed to be diagonal. In this case the two
correlations can be calibrated separately and it
is refered to as diag-diag calibration.
The second case is when the Jones matrices are
assumed to be diagonal but the visibility data
are full \(2\times 2\) matrices. This is
refered to as diag calibration. The final
scenario is when both the full \(2\times 2\)
Jones matrices and the full \(2\times 2\)
visibilities are used for calibration. This is
simply refered to as calibration. The specific
scenario is determined from the shapes of the input
gains and the input data.
This module also provides a number of utilities which are useful for calibration.
Utils¶
Numpy¶
corrupt_vis(time_bin_indices, …) |
Corrupts model visibilities with arbitrary Jones terms. |
residual_vis(time_bin_indices, …) |
Computes residual visibilities given model visibilities and gains solutions. |
correct_vis(time_bin_indices, …) |
Apply inverse of direction independent gains to visibilities to generate corrected visibilities. |
compute_and_corrupt_vis(time_bin_indices, …) |
Corrupts time variable component model with arbitrary Jones terms. |
-
africanus.calibration.utils.corrupt_vis(time_bin_indices, time_bin_counts, antenna1, antenna2, jones, model)[source]¶ Corrupts model visibilities with arbitrary Jones terms.
Parameters: time_bin_indices :
numpy.ndarrayThe start indices of the time bins of shape
(utime)time_bin_counts :
numpy.ndarrayThe counts of unique time in each time bin of shape
(utime)antenna1 :
numpy.ndarrayFirst antenna indices of shape
(row,).antenna2 :
numpy.ndarraySecond antenna indices of shape
(row,)jones :
numpy.ndarrayGains of shape
(time, ant, chan, dir, corr)or(time, ant, chan, dir, corr, corr).model :
numpy.ndarrayModel data values of shape
(row, chan, dir, corr)or(row, chan, dir, corr, corr).Returns: vis :
numpy.ndarrayvisibilities of shape
(time, ant, chan, dir, corr)or(time, ant, chan, dir, corr, corr).
-
africanus.calibration.utils.residual_vis(time_bin_indices, time_bin_counts, antenna1, antenna2, jones, vis, flag, model)[source]¶ Computes residual visibilities given model visibilities and gains solutions.
Parameters: time_bin_indices :
numpy.ndarrayThe start indices of the time bins of shape
(utime)time_bin_counts :
numpy.ndarrayThe counts of unique time in each time bin of shape
(utime)antenna1 :
numpy.ndarrayFirst antenna indices of shape
(row,).antenna2 :
numpy.ndarraySecond antenna indices of shape
(row,)jones :
numpy.ndarrayGain solutions of shape
(time, ant, chan, dir, corr)or(time, ant, chan, dir, corr, corr).vis :
numpy.ndarrayData values of shape
(row, chan, corr). or(row, chan, corr, corr).flag :
numpy.ndarrayFlag data of shape
(row, chan, corr)or(row, chan, corr, corr)model :
numpy.ndarrayModel data values of shape
(row, chan, dir, corr)or(row, chan, dir, corr, corr).Returns: residual :
numpy.ndarrayResidual visibilities of shape
(time, ant, chan, dir, corr)or(time, ant, chan, dir, corr, corr).
-
africanus.calibration.utils.correct_vis(time_bin_indices, time_bin_counts, antenna1, antenna2, jones, vis, flag)[source]¶ Apply inverse of direction independent gains to visibilities to generate corrected visibilities. For a measurement model of the form
\[V_{pq} = G_{p} X_{pq} G_{q}^H + n_{pq}\]the corrected visibilities are defined as
\[C_{pq} = G_{p}^{-1} V_{pq} G_{q}^{-H}\]The corrected visibilities therefore have a non-trivial noise contribution. Note it is only possible to form corrected data from direction independent gains solutions so the
diraxis on the jones terms should always be one.Parameters: time_bin_indices :
numpy.ndarrayThe start indices of the time bins of shape
(utime).time_bin_counts :
numpy.ndarrayThe counts of unique time in each time bin of shape
(utime).antenna1 :
numpy.ndarrayAntenna 1 index used to look up the antenna Jones for a particular baseline with shape
(row,).antenna2 :
numpy.ndarrayAntenna 2 index used to look up the antenna Jones for a particular baseline with shape
(row,).jones :
numpy.ndarrayGain solutions of shape
(time, ant, chan, dir, corr)or(time, ant, chan, dir, corr, corr).vis :
numpy.ndarrayData values of shape
(row, chan, corr)or(row, chan, corr, corr).flag :
numpy.ndarrayFlag data of shape
(row, chan, corr)or(row, chan, corr, corr).Returns
——-
corrected_vis :
numpy.ndarrayTrue visibilities of shape
(row,chan,corr_1,corr_2)
-
africanus.calibration.utils.compute_and_corrupt_vis(time_bin_indices, time_bin_counts, antenna1, antenna2, jones, model, uvw, freq, lm)[source]¶ Corrupts time variable component model with arbitrary Jones terms. Currrently only time variable point source models are supported.
Parameters: time_bin_indices :
numpy.ndarrayThe start indices of the time bins of shape
(utime)time_bin_counts :
numpy.ndarrayThe counts of unique time in each time bin of shape
(utime)antenna1 :
numpy.ndarrayFirst antenna indices of shape
(row,).antenna2 :
numpy.ndarraySecond antenna indices of shape
(row,)jones :
numpy.ndarrayGains of shape
(utime, ant, chan, dir, corr)or(utime, ant, chan, dir, corr, corr).model :
numpy.ndarrayModel image as a function of time with shape
(utime, chan, dir, corr)or(utime, chan, dir, corr, corr).uvw :
numpy.ndarrayuvw coordinates of shape
(row, 3)lm :
numpy.ndarraySource lm coordinates as a function of time
(utime, dir, 2)Returns: vis :
numpy.ndarrayvisibilities of shape
(row, chan, corr)or(row, chan, corr, corr).
Dask¶
corrupt_vis(time_bin_indices, …) |
Corrupts model visibilities with arbitrary Jones terms. |
residual_vis(time_bin_indices, …) |
Computes residual visibilities given model visibilities and gains solutions. |
correct_vis(time_bin_indices, …) |
Apply inverse of direction independent gains to visibilities to generate corrected visibilities. |
compute_and_corrupt_vis(time_bin_indices, …) |
Corrupts time variable component model with arbitrary Jones terms. |
-
africanus.calibration.utils.dask.corrupt_vis(time_bin_indices, time_bin_counts, antenna1, antenna2, jones, model)[source]¶ Corrupts model visibilities with arbitrary Jones terms.
Parameters: time_bin_indices :
dask.array.ArrayThe start indices of the time bins of shape
(utime)time_bin_counts :
dask.array.ArrayThe counts of unique time in each time bin of shape
(utime)antenna1 :
dask.array.ArrayFirst antenna indices of shape
(row,).antenna2 :
dask.array.ArraySecond antenna indices of shape
(row,)jones :
dask.array.ArrayGains of shape
(time, ant, chan, dir, corr)or(time, ant, chan, dir, corr, corr).model :
dask.array.ArrayModel data values of shape
(row, chan, dir, corr)or(row, chan, dir, corr, corr).Returns: vis :
dask.array.Arrayvisibilities of shape
(time, ant, chan, dir, corr)or(time, ant, chan, dir, corr, corr).
-
africanus.calibration.utils.dask.residual_vis(time_bin_indices, time_bin_counts, antenna1, antenna2, jones, vis, flag, model)[source]¶ Computes residual visibilities given model visibilities and gains solutions.
Parameters: time_bin_indices :
dask.array.ArrayThe start indices of the time bins of shape
(utime)time_bin_counts :
dask.array.ArrayThe counts of unique time in each time bin of shape
(utime)antenna1 :
dask.array.ArrayFirst antenna indices of shape
(row,).antenna2 :
dask.array.ArraySecond antenna indices of shape
(row,)jones :
dask.array.ArrayGain solutions of shape
(time, ant, chan, dir, corr)or(time, ant, chan, dir, corr, corr).vis :
dask.array.ArrayData values of shape
(row, chan, corr). or(row, chan, corr, corr).flag :
dask.array.ArrayFlag data of shape
(row, chan, corr)or(row, chan, corr, corr)model :
dask.array.ArrayModel data values of shape
(row, chan, dir, corr)or(row, chan, dir, corr, corr).Returns: residual :
dask.array.ArrayResidual visibilities of shape
(time, ant, chan, dir, corr)or(time, ant, chan, dir, corr, corr).
-
africanus.calibration.utils.dask.correct_vis(time_bin_indices, time_bin_counts, antenna1, antenna2, jones, vis, flag)[source]¶ Apply inverse of direction independent gains to visibilities to generate corrected visibilities. For a measurement model of the form
\[V_{pq} = G_{p} X_{pq} G_{q}^H + n_{pq}\]the corrected visibilities are defined as
\[C_{pq} = G_{p}^{-1} V_{pq} G_{q}^{-H}\]The corrected visibilities therefore have a non-trivial noise contribution. Note it is only possible to form corrected data from direction independent gains solutions so the
diraxis on the jones terms should always be one.Parameters: time_bin_indices :
dask.array.ArrayThe start indices of the time bins of shape
(utime).time_bin_counts :
dask.array.ArrayThe counts of unique time in each time bin of shape
(utime).antenna1 :
dask.array.ArrayAntenna 1 index used to look up the antenna Jones for a particular baseline with shape
(row,).antenna2 :
dask.array.ArrayAntenna 2 index used to look up the antenna Jones for a particular baseline with shape
(row,).jones :
dask.array.ArrayGain solutions of shape
(time, ant, chan, dir, corr)or(time, ant, chan, dir, corr, corr).vis :
dask.array.ArrayData values of shape
(row, chan, corr)or(row, chan, corr, corr).flag :
dask.array.ArrayFlag data of shape
(row, chan, corr)or(row, chan, corr, corr).Returns
——-
corrected_vis :
dask.array.ArrayTrue visibilities of shape
(row,chan,corr_1,corr_2)
-
africanus.calibration.utils.dask.compute_and_corrupt_vis(time_bin_indices, time_bin_counts, antenna1, antenna2, jones, model, uvw, freq, lm)[source]¶ Corrupts time variable component model with arbitrary Jones terms. Currrently only time variable point source models are supported.
Parameters: time_bin_indices :
dask.array.ArrayThe start indices of the time bins of shape
(utime)time_bin_counts :
dask.array.ArrayThe counts of unique time in each time bin of shape
(utime)antenna1 :
dask.array.ArrayFirst antenna indices of shape
(row,).antenna2 :
dask.array.ArraySecond antenna indices of shape
(row,)jones :
dask.array.ArrayGains of shape
(utime, ant, chan, dir, corr)or(utime, ant, chan, dir, corr, corr).model :
dask.array.ArrayModel image as a function of time with shape
(utime, chan, dir, corr)or(utime, chan, dir, corr, corr).uvw :
dask.array.Arrayuvw coordinates of shape
(row, 3)lm :
dask.array.ArraySource lm coordinates as a function of time
(utime, dir, 2)Returns: vis :
dask.array.Arrayvisibilities of shape
(row, chan, corr)or(row, chan, corr, corr).
Phase only¶
Numpy¶
compute_jhr(time_bin_indices, …) |
Computes the residual projected in to gain space. |
compute_jhj(time_bin_indices, …) |
Computes the diagonal of the Hessian required to perform phase-only maximum likelihood calibration. |
compute_jhj_and_jhr(time_bin_indices, …) |
Computes the diagonal of the Hessian and the residual locally projected in to gain space. |
gauss_newton(time_bin_indices, …[, tol, …]) |
Performs phase-only maximum likelihood calibration using a Gauss-Newton optimisation algorithm. |
-
africanus.calibration.phase_only.compute_jhr(time_bin_indices, time_bin_counts, antenna1, antenna2, jones, residual, model, flag)[source]¶ Computes the residual projected in to gain space.
Parameters: time_bin_indices :
numpy.ndarrayThe start indices of the time bins of shape
(utime)time_bin_counts :
numpy.ndarrayThe counts of unique time in each time bin of shape
(utime)antenna1 :
numpy.ndarrayFirst antenna indices of shape
(row,).antenna2 :
numpy.ndarraySecond antenna indices of shape
(row,)jones :
numpy.ndarrayGain solutions of shape
(time, ant, chan, dir, corr)or(time, ant, chan, dir, corr, corr).residual :
numpy.ndarrayResidual values of shape
(row, chan, corr). or(row, chan, corr, corr).model :
numpy.ndarrayModel data values of shape
(row, chan, dir, corr)or(row, chan, dir, corr, corr).flag :
numpy.ndarrayFlag data of shape
(row, chan, corr)or(row, chan, corr, corr)Returns: jhr :
numpy.ndarrayThe residual projected into gain space shape
(time, ant, chan, dir, corr)or(time, ant, chan, dir, corr, corr).
-
africanus.calibration.phase_only.compute_jhj(time_bin_indices, time_bin_counts, antenna1, antenna2, jones, model, flag)[source]¶ Computes the diagonal of the Hessian required to perform phase-only maximum likelihood calibration. Currently assumes scalar or diagonal inputs.
Parameters: time_bin_indices :
numpy.ndarrayThe start indices of the time bins of shape
(utime)time_bin_counts :
numpy.ndarrayThe counts of unique time in each time bin of shape
(utime)antenna1 :
numpy.ndarrayFirst antenna indices of shape
(row,).antenna2 :
numpy.ndarraySecond antenna indices of shape
(row,)jones :
numpy.ndarrayGain solutions of shape
(time, ant, chan, dir, corr)or(time, ant, chan, dir, corr, corr).model :
numpy.ndarrayModel data values of shape
(row, chan, dir, corr)or(row, chan, dir, corr, corr).flag :
numpy.ndarrayFlag data of shape
(row, chan, corr)or(row, chan, corr, corr)Returns: jhj :
numpy.ndarrayThe diagonal of the Hessian of shape
(time, ant, chan, dir, corr)or(time, ant, chan, dir, corr, corr).
-
africanus.calibration.phase_only.compute_jhj_and_jhr(time_bin_indices, time_bin_counts, antenna1, antenna2, jones, residual, model, flag)[source]¶ Computes the diagonal of the Hessian and the residual locally projected in to gain space.
Parameters: time_bin_indices :
numpy.ndarrayThe start indices of the time bins of shape
(utime)time_bin_counts :
numpy.ndarrayThe counts of unique time in each time bin of shape
(utime)antenna1 :
numpy.ndarrayFirst antenna indices of shape
(row,).antenna2 :
numpy.ndarraySecond antenna indices of shape
(row,)jones :
numpy.ndarrayGain solutions of shape
(time, ant, chan, dir, corr)or(time, ant, chan, dir, corr, corr).residual :
numpy.ndarrayResidual values of shape
(row, chan, corr). or(row, chan, corr, corr).model :
numpy.ndarrayModel data values of shape
(row, chan, dir, corr)or(row, chan, dir, corr, corr).flag :
numpy.ndarrayFlag data of shape
(row, chan, corr)or(row, chan, corr, corr)Returns: jhj :
numpy.ndarrayThe diagonal of the Hessian of shape
(time, ant, chan, dir, corr)or(time, ant, chan, dir, corr, corr).jhr :
numpy.ndarrayResiduals projected into signal space of shape
(time, ant, chan, dir, corr)or(time, ant, chan, dir, corr, corr).
-
africanus.calibration.phase_only.gauss_newton(time_bin_indices, time_bin_counts, antenna1, antenna2, jones, vis, flag, model, weight, tol=0.0001, maxiter=100)[source]¶ Performs phase-only maximum likelihood calibration using a Gauss-Newton optimisation algorithm. Currently only DIAG mode is supported.
Parameters: time_bin_indices :
numpy.ndarrayThe start indices of the time bins of shape
(utime)time_bin_counts :
numpy.ndarrayThe counts of unique time in each time bin of shape
(utime)antenna1 :
numpy.ndarrayFirst antenna indices of shape
(row,).antenna2 :
numpy.ndarraySecond antenna indices of shape
(row,).jones :
numpy.ndarrayGain solutions of shape
(time, ant, chan, dir, corr)or(time, ant, chan, dir, corr, corr).vis :
numpy.ndarrayData values of shape
(row, chan, corr)or(row, chan, corr, corr).flag :
numpy.ndarrayFlag data of shape
(row, chan, corr)or(row, chan, corr, corr).model :
numpy.ndarrayModel data values of shape
(row, chan, dir, corr)or(row, chan, dir, corr, corr).weight :
numpy.ndarrayWeight spectrum of shape
(row, chan, corr). If the channel axis is missing weights are duplicated for each channel.tol: float, optional
The tolerance of the solver. Defaults to 1e-4.
maxiter: int, optional
The maximum number of iterations. Defaults to 100.
Returns: gains :
numpy.ndarrayGain solutions of shape
(time, ant, chan, dir, corr)or shape(time, ant, chan, dir, corr, corr)jhj :
numpy.ndarrayThe diagonal of the Hessian of shape
(time, ant, chan, dir, corr)or shape(time, ant, chan, dir, corr, corr)jhr :
numpy.ndarrayResiduals projected into gain space of shape
(time, ant, chan, dir, corr)or shape(time, ant, chan, dir, corr, corr).k: int
Number of iterations (will equal maxiter if not converged)
Dask¶
compute_jhr(time_bin_indices, …) |
Computes the residual projected in to gain space. |
compute_jhj(time_bin_indices, …) |
Computes the diagonal of the Hessian required to perform phase-only maximum likelihood calibration. |
-
africanus.calibration.phase_only.dask.compute_jhr(time_bin_indices, time_bin_counts, antenna1, antenna2, jones, residual, model, flag)[source]¶ Computes the residual projected in to gain space.
Parameters: time_bin_indices :
dask.array.ArrayThe start indices of the time bins of shape
(utime)time_bin_counts :
dask.array.ArrayThe counts of unique time in each time bin of shape
(utime)antenna1 :
dask.array.ArrayFirst antenna indices of shape
(row,).antenna2 :
dask.array.ArraySecond antenna indices of shape
(row,)jones :
dask.array.ArrayGain solutions of shape
(time, ant, chan, dir, corr)or(time, ant, chan, dir, corr, corr).residual :
dask.array.ArrayResidual values of shape
(row, chan, corr). or(row, chan, corr, corr).model :
dask.array.ArrayModel data values of shape
(row, chan, dir, corr)or(row, chan, dir, corr, corr).flag :
dask.array.ArrayFlag data of shape
(row, chan, corr)or(row, chan, corr, corr)Returns: jhr :
dask.array.ArrayThe residual projected into gain space shape
(time, ant, chan, dir, corr)or(time, ant, chan, dir, corr, corr).
-
africanus.calibration.phase_only.dask.compute_jhj(time_bin_indices, time_bin_counts, antenna1, antenna2, jones, model, flag)[source]¶ Computes the diagonal of the Hessian required to perform phase-only maximum likelihood calibration. Currently assumes scalar or diagonal inputs.
Parameters: time_bin_indices :
dask.array.ArrayThe start indices of the time bins of shape
(utime)time_bin_counts :
dask.array.ArrayThe counts of unique time in each time bin of shape
(utime)antenna1 :
dask.array.ArrayFirst antenna indices of shape
(row,).antenna2 :
dask.array.ArraySecond antenna indices of shape
(row,)jones :
dask.array.ArrayGain solutions of shape
(time, ant, chan, dir, corr)or(time, ant, chan, dir, corr, corr).model :
dask.array.ArrayModel data values of shape
(row, chan, dir, corr)or(row, chan, dir, corr, corr).flag :
dask.array.ArrayFlag data of shape
(row, chan, corr)or(row, chan, corr, corr)Returns: jhj :
dask.array.ArrayThe diagonal of the Hessian of shape
(time, ant, chan, dir, corr)or(time, ant, chan, dir, corr, corr).