nenupy.instru.interferometer


Interferometric Array

class nenupy.instru.interferometer.Baseline(itrf_positions)[source]

Bases: object

Class to handle interferometric baseline operations.

property distance
property flatten
class nenupy.instru.interferometer.Interferometer(position, antenna_names, antenna_positions, antenna_gains, antenna_delays=None, antenna_weights=None)[source]

Bases: ABC

Abstract base class for all phased-array/interferometer classes.

Attributes Summary

antenna_gains

Antenna gains.

antenna_names

Antenna names.

antenna_positions

Antenna positions.

baselines

Instrument baselines.

position

Array's position.

size

Number of elements belonging to the array.

antenna_weights

antenna_delays

Add delay errors for each antennae They could be cable connection errors during construction, cables of wrong length, ...

Methods Summary

angular_resolution([frequency])

Computes the angular resolution of the antenna array.

array_factor(sky, pointing[, ...])

Computes the array factor of the antenna distribution.

beam(sky, pointing[, return_complex, normalize])

Computes the phased-array response \(\mathcal{G}\) over the sky for a given pointing.

confusion_noise([frequency, lofar])

Confusion rms noise \(\sigma_{\rm c}\) (parameter used for specifying the width of the confusion distribution) computed as:

plot(**kwargs)

Plots the antenna distribution.

sefd([frequency, elevation, efficiency, ...])

Computes the System Equivalent Flux Density (SEFD or system sensitivity).

sensitivity([frequency, mode, dt, df, ...])

Computes the sensititivy of the array with respect to the observing configuration.

system_temperature([frequency, ...])

Computes the System Noise Temperature \(T_{\rm sys}\).

Attributes and Methods Documentation

angular_resolution(frequency=<Quantity 50. MHz>)[source]

Computes the angular resolution of the antenna array.

The full width at half maximum (FWHM) \(\theta\) is approximated as follows:

\[\theta = \frac{\lambda}{D}\]

where \(\lambda\) is the wavelength and \(D\) is is the length of the maximum physical separation of the antennas in the array.

Parameters:

frequency (Quantity) – Frequency at which the angular resolution is evaluated.

Returns:

Angular resolution (FWHM) of the instrument.

Return type:

Quantity

Example:
>>> import astropy.units as u
>>> <instrument>.angular_resolution(
        frequency=50*u.MHz
    )
property antenna_delays

Add delay errors for each antennae They could be cable connection errors during construction, cables of wrong length, …

property antenna_gains

Antenna gains. This is an array of callable (methods or functions) defining the radiation pattern of each antenna.

Setter:

Array of antenna gains.

Getter:

Array of antenna gains.

Type:

ndarray of callable

property antenna_names

Antenna names.

Setter:

Array of antenna names.

Getter:

Array of antenna names.

Type:

ndarray

property antenna_positions

Antenna positions. The positions should be shaped as (n_ant, 3)

Setter:

Array of antenna positions.

Getter:

Array of antenna positions.

Type:

ndarray

property antenna_weights
array_factor(sky, pointing, return_complex=False, normalize=True)[source]

Computes the array factor of the antenna distribution.

\[\mathcal{F}(\nu, \phi, \theta) = \sum_{\rm ant} w_{\rm ant} e^{ i \mathbf{k}(\nu, \phi, \theta) \cdot \mathbf{r}_{\rm ant}}\]

where \(\mathbf{k} = \frac{2\pi}{\lambda} (\cos \phi \cos \theta, \sin \phi \cos \theta, \sin \theta )\) is the wave vector for a wave propagation in a direction described by spherical coordinates, \(\lambda\) is the wavelength, \(\phi\) is the azimuth, \(\theta\) is the elevation, \(\mathbf{r}_{\rm ant}\) is the antenna position matrix, \(w_{\rm ant}\) is the weight of the antenna (defined in feed_weights).

This method considers the sky as the desired output (in terms of time, frequency and sky positions). It evaluates the effective pointing directions for every time step defined in sky regarding the pointing input.

Parameters:
  • sky (Sky) – Desired output contained in a Sky instance. (time, frequency and coordinates are used as inputs for the computation).

  • pointing (Pointing) – Instance of Pointing that defines the targeted pointing directions over the time.

  • return_complex (bool) – Return complex array factor if True or power if False

  • normalize (bool) – Return the normalized array factor. Default is True.

Returns:

Array factor of the antenna distribution shaped as (time, frequency, 1, coordinates).

Return type:

Array

See also

Sky and Pointing

property baselines

Instrument baselines.

Getter:

Baselines.

Type:

Baseline

beam(sky, pointing, return_complex=False, normalize=True)[source]

Computes the phased-array response \(\mathcal{G}\) over the sky for a given pointing.

\[\mathcal{G}(\nu, \phi, \theta) = \sum_{\rm ant} \mathcal{F}(\nu, \phi, \theta) \mathcal{G}_{\rm ant} (\nu, \phi, \theta)\]

where \(\nu\) is the frequency, \(\phi\) is the azimuth, \(\theta\) is the elevation, \(\mathcal{G}_{\rm ant}\) is the individual array element radiation pattern and \(\mathcal{F}\) is the array factor.

This method considers the sky as the desired output (in terms of time, frequency, polarization and sky positions). It evaluates the effective pointing directions for every time step defined in sky regarding the pointing input.

Parameters:
  • sky (Sky) – Desired output contained in a Sky instance.

  • pointing (Pointing) – Instance of Pointing that defines the targeted pointing directions over the time.

Returns:

The instance of Sky given as input is returned, its attribute value is updated with the result of the beam computation (stored as an Array) and shaped as (time, frequency, polarization, coordinates).

Return type:

Sky

See also

array_factor()

confusion_noise(frequency=<Quantity 50. MHz>, lofar=True)[source]

Confusion rms noise \(\sigma_{\rm c}\) (parameter used for specifying the width of the confusion distribution) computed as:

\[\left( \frac{\sigma_{\rm c}}{\rm{mJy}\, \rm{beam}^{-1}} \right) \simeq 0.2 \left( \frac{\nu}{\rm GHz} \right)^{-0.7} \left( \frac{\theta}{\rm arcmin} \right)^{2}\]

or (if lofar=True):

\[\left( \frac{\sigma_{\rm c}}{\mu\rm{Jy}\, \rm{beam}^{-1}} \right) \simeq 30 \left( \frac{\nu}{74 {\rm MHz}} \right)^{-0.7} \left( \frac{\theta}{\rm arcsec} \right)^{1.54}\]

where \(\nu\) is the frequency and \(\theta\) is the radiotelescope FWHM (see angular_resolution()).

Individual sources fainter than about \(5\sigma_{\rm c}\) cannot be detected reliably.

Parameters:
  • freq (float or Quantity) – Frequency at which computing the confusion noise. In MHz if no unit is provided. Default is 50 MHz.

  • miniarrays (int, list or ndarray) – Mini-Array indices to take into account. Default is None (all available MAs).

  • lofar – If set to True (recommended), the confusion noise is estimated using Eq. 6 of van Haarlem et al. (2013).

Type:

bool

Returns:

Confusion rms noise in Jy/beam

Return type:

Quantity

Example:
>>> import astropy.units as u
>>> <instrument>.confusion_noise(
        frequency=50*u.MHz
    )
abstract effective_area(frequency=<Quantity 50. MHz>, elevation=<Quantity 90. deg>)[source]

This method needs to be defined in child classes.

abstract instrument_temperature(frequency=<Quantity 50. MHz>, **kwargs)[source]

This method needs to be defined in child classes.

plot(**kwargs)[source]

Plots the antenna distribution.

Parameters:
  • figsize (tuple) – Size of the figure. Default is (10, 10).

  • figname (str) – File name of the figure to save. Default is '', i.e. show the figure without saving it.

  • xlim (tuple) – X-axis limits. Default is auto-scaling.

  • ylim (tuple) – Y-axis limits. Default is auto-scaling.

  • show_names (bool) – Print the antenna names. Default is True.

  • patches (tuple`(`list, colors) of Polygons) – Matplotlib Polygons

property position

Array’s position.

Setter:

Position of the array.

Getter:

Position of the array.

Type:

EarthLocation

sefd(frequency=<Quantity 50. MHz>, elevation=<Quantity 90. deg>, efficiency=1.0, decoherence=1.0, source_spectrum={}, **kwargs)[source]

Computes the System Equivalent Flux Density (SEFD or system sensitivity).

\[S_{\rm sys} = \xi \frac{2 k_{\rm B}}{ \eta A_{\rm eff}(\nu, \theta)} T_{\rm sys} (\nu)\]

with \(T_{\rm sys}\) the system_temperature(), the efficiency \(\eta\), \(\nu\) the frequency, \(\theta\) the elevation, \(\xi\) the decoherence factor, and \(k_{\rm B}\) the Boltzmann constant.

Parameters:
  • frequency (Quantity) – Frequency at which the SEFD will be computed. If an array is given as input, the output will be of same shape. Default if 50 MHz.

  • elevation (Quantity) – Pointing elevation impacting the effective_area(). Default is 90 deg.

  • efficiency (float) – Effective area reducing factor. Default is 1., it cannot be greater than 1..

  • decoherence (float) – Parameter that reflects other uncertainties (particularly the unperfect phasing system). Default is 1..

  • source_spectrum (dict of callable) – By default the system temperature is evaluated using a mean Galactic temperature. However, if a bright source is targeted, the noise introduced can be under-estimated. Therefore, one can provide a callable object that takes as inputs a frequency array (of type Quantity) and returns the source flux density in Jansky (of type (as Quantity).

Returns:

SEFD in Janskys.

Return type:

Quantity

sensitivity(frequency=<Quantity 50. MHz>, mode=ObservingMode.BEAMFORMING, dt=<Quantity 1. s>, df=<Quantity 195.3125 kHz>, elevation=<Quantity 90. deg>, efficiency=1.0, decoherence=1.0, source_spectrum={}, **kwargs)[source]

Computes the sensititivy of the array with respect to the observing configuration. The sensitivity computation depends on the observing mode of the instrument:

  • for the imaging mode:

    \[\sigma_{\rm im} = \frac{S_{\rm sys}(\nu, \theta, \eta, \xi)}{ \sqrt{N(N-1) 2 \Delta \nu\, \Delta t} }\]
  • for the beamforming mode:

    \[\sigma_{\rm bf} = \frac{S_{\rm sys}(\nu, \theta, \eta, \xi)}{ \sqrt{2 \Delta \nu\, \Delta t} }\]

where \(\nu\) is the frequency, \(\theta\) is the elevation, \(\eta\) is the effective area efficiency, \(\xi\) is the decoherence factor, \(\Delta t\) is the integration time, \(\Delta \nu\) is the bandwidth, \(N\) is the antenna number, and \(S_{\rm sys}\) is the System Equivalent Flux Density (which also depends on the source_spectrum argument, see sefd()).

Parameters:
  • frequency (Quantity) – Frequency at which the sensitivity will be evaluated. If an array is given as input, the output will be of same shape. Default if 50 MHz.

  • mode (ObservingMode) – Observing mode, either ObservingMode.BEAMFORMING or ObservingMode.IMAGING, default is the former.

  • dt (Quantity) – Integration time. Default is 1 sec.

  • df (Quantity) – Observing bandwidth. Default is 195.3125 kHz.

  • elevation (Quantity) – Pointing elevation impacting the effective_area(). Default is 90 deg.

  • efficiency (float) – Effective area reducing factor. Default is 1., it cannot be greater than 1..

  • decoherence (float) – Parameter that reflects other uncertainties (particularly the unperfect phasing system). Default is 1..

  • source_spectrum (dict of callable) – By default the system temperature is evaluated using a mean Galactic temperature. However, if a bright source is targeted, the noise introduced can be under-estimated. Therefore, one can provide a callable object that takes as inputs a frequency array (of type Quantity) and returns the source flux density in Jansky (of type Quantity).

Returns:

Array sensitivity.

Return type:

Quantity

Example:
>>> from nenupy.instru.interferometer import ObservingMode
>>> import astropy.units as u
>>> <instrument>.sensitivity(
        frequency=50*u.MHz,
        mode=ObservingMode.IMAGING,
        dt=1*u.s,
        df=3*u.kHz
    )
property size

Number of elements belonging to the array.

Getter:

Size of the array.

Type:

int

system_temperature(frequency=<Quantity 50. MHz>, source_spectrum={}, efficiency=1.0, elevation=<Quantity 90. deg>, **kwargs)[source]

Computes the System Noise Temperature \(T_{\rm sys}\). It is computed as follows:

\[T_{\rm sys} = T_{\rm sky} + T_{\rm inst} + \sum_{\rm src} T_{\rm src}\]

where \(T_{\rm sky}\) is an approximation of the low-frequency sky temperature dominated by Galactic emission and \(T_{\rm inst}\) is the instrumental noise temperature (which depends on the current instrument instance). \(T_{\rm src}\) is the antenna temperature induced by a given source whose spectrum is defined in the source_spectrum argument computed as:

\[T_{\rm src} = \frac{F_{\rm src} \eta A_{\rm eff}}{2 k_{\rm B}}\]

where \(F_{\rm src}\) is the source spectrum, \(\eta\) is the efficiency of the effective area \(A_{\rm eff}\).

Parameters:
  • frequency (Quantity) – Frequency for the System Temperature computation. Default is 50 MHz.

  • elevation (Quantity) – Pointing elevation impacting the effective_area(). Default is 90 deg.

  • efficiency (float) – Effective area reducing factor. Default is 1., it cannot be greater than 1..

  • source_spectrum (dict of callable) – By default the system temperature is evaluated using a mean Galactic temperature. However, if a bright source is targeted, the noise introduced can be under-estimated. Therefore, one can provide a callable object that takes as inputs a frequency array (of type Quantity) and returns the source flux density in Jansky (of type Quantity).

Returns:

System Temperature in Kelvins.

Return type:

Quantity

class nenupy.instru.interferometer.ObservingMode(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)[source]

Bases: Enum

Enumerator of the different available observing modes of NenuFAR.

BEAMFORMING = 1
IMAGING = 2