nenupy.instru.nenufar


NenuFAR Array Classes

Inheritance diagram of nenupy.instru.nenufar.MiniArray, nenupy.instru.nenufar.NenuFAR

MiniArray([index, antenna_delays, ...])

Main class to handle a NenuFAR Mini-Array antenna distribution.

NenuFAR([miniarray_antennas, include_remote_mas])

Main class to handle a NenuFAR array.

class nenupy.instru.nenufar.MiniArray(index=0, antenna_delays=None, antenna_weights=None)[source]

Bases: Interferometer

Main class to handle a NenuFAR Mini-Array antenna distribution.

Added in version 2.0.0.

Parameters:

index (int) – Mini-Array index. ‘Core’ Mini-Arrays have indices ranging from 0 to 95. ‘Remote’ Mini-Arrays have indices ranging from 100 to 105.

Example:

Instantiating MiniArray:

>>> from nenupy.instru import MiniArray
>>> ma = MiniArray(index=0)

Sub-arraying on an existing MiniArray instance:

>>> sub_ma = ma["Ant01", "Ant06", "Ant11"]
>>> sub_ma.antenna_names
array(['Ant01', 'Ant06', 'Ant11'], dtype='<U5')

Using slice object (converted in ndarray using r_):

>>> import numpy as np
>>> sub_ma = ma[np.r_[2:10]]
>>> sub_ma.size
8

Combining two MiniArray instances:

>>> ma1 = MiniArray(index=0)["Ant01", "Ant06"]
>>> ma2 = MiniArray(index=0)["Ant08", "Ant12"]
>>> combined_ma = ma1 + ma2
>>> combined_ma.antenna_names
array(['Ant01', 'Ant06', 'Ant08', 'Ant12'], dtype='<U5')

See also

More details on this class usage can be found in Array Configuration and Instrument Properties.

Attributes Summary

index

Mini-Array index.

rotation

Mini-Array rotation.

position

Array's position.

antenna_names

Antenna names.

antenna_positions

Antenna positions.

antenna_gains

Antenna gains.

baselines

Instrument baselines.

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

beam(sky, pointing[, configuration, ...])

Computes the Mini-Array beam over the sky for a given pointing.

effective_area([frequency, elevation])

Computes the effective area of a NenuFAR Mini-Array.

instrument_temperature([frequency, lna_filter])

Instrument temperature at a given frequency.

attenuation_from_zenith(coordinates[, time, ...])

Returns the attenuation factor evaluated at given coordinates compared to the zenithal Mini-Array beam gain.

analog_pointing(pointing, configuration)

Converts the desired pointing to the effective pointing which depends on the available pointing positions defined on a grid due to analog cable delays.

beamsquint_correction(coords[, frequency])

Corrects for the beamsquint effect.

plot(**kwargs)

Plots the antenna distribution.

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

Computes the array factor of the antenna distribution.

system_temperature([frequency, ...])

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

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.

angular_resolution([frequency])

Computes the angular resolution of the antenna array.

confusion_noise([frequency, lofar])

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

Attributes and Methods Documentation

analog_pointing(pointing, configuration)[source]

Converts the desired pointing to the effective pointing which depends on the available pointing positions defined on a grid due to analog cable delays.

attenuation_from_zenith(coordinates, time=<Time object: scale='utc' format='datetime' value=2024-05-03 08:29:42.563431>, frequency=<Quantity 50. MHz>, polarization=Polarization.NW)[source]

Returns the attenuation factor evaluated at given coordinates compared to the zenithal Mini-Array beam gain.

Parameters:
  • coordinates (SkyCoord) – Sky positions equatorial coordinates.

  • time (Time) – UTC time at which the attenuation is evaluated. Default is now.

  • frequency – Frequency at which the attenuation is evaluated. Default is 50 MHz.

  • polarization (Polarization) – NenuFAR antenna polarization. Default is Polarization.NW.

Returns:

Attenuation factor shaped as (time, frequency, polarization, coordinates). NaN is returned for any coordinates that is below the horizon.

Return type:

ndarray

Example:
>>> from nenupy.instru.nenufar import MiniArray
>>> from astropy.coordinates import SkyCoord
>>> ma = MiniArray(index=0)
>>> attenuation = ma.attenuation_from_zenith(
        coordinates=SkyCoord.from_name("Cyg A")
    )
>>> from nenupy.instru.nenufar import MiniArray
>>> from astropy.coordinates import SkyCoord
>>> import astropy.units as u
>>> ma = MiniArray(index=0)
>>> attenuation = ma.attenuation_from_zenith(
        coordinates=SkyCoord.from_name("Cyg A"),
        frequency=np.linspace(20, 80, 10)*u.MHz
    )

Added in version 2.0.0.

beam(sky, pointing, configuration=<nenupy.instru.nenufar.NenuFAR_Configuration object>, return_complex=False, normalize=True)[source]

Computes the Mini-Array beam over the sky for a given pointing.

\[\mathcal{G}_{\rm MA}(\nu, \phi, \theta) = \mathcal{F}_{\rm MA}(\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 NenuFAR dipole antenna radiation pattern and \(\mathcal{F}_{\rm MA}\) 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. (time, frequency, polarization and coordinates are used as inputs for the computation).

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

  • configuration (NenuFAR_Configuration) – NenuFAR configuration to consider during the beam simulation. The beamsquint correction and its frequency setting are defined here. Default is NenuFAR_Configuration(beamsquint_correction=True, beamsquint_frequency=50MHz).

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

Example:

Load the required librairies:

>>> from nenupy.instru import MiniArray, Polarization
>>> from nenupy.astro.sky import HpxSky
>>> from nenupy.astro.pointing import Pointing
>>> import astropy.units as u
>>> from astropy.time import Time, TimeDelta

Define a desired Sky output:

>>> sky = HpxSky(
        resolution=1.*u.deg,
        frequency=np.array([25, 50, 75])*u.MHz,
        polarization=np.array([Polarization.NW, Polarization.NE]),
        time=Time("2021-10-15 20:00:00")
    )

Define the pointing of the Mini-Array:

>>> ma_pointing = Pointing.zenith_tracking(
        time=Time("2021-10-15 00:00:00"),
        duration=TimeDelta(3600*24, format="sec")
    )

Select the Mini-Array (and possibly its antenna distribution) and compute its response pattern:

>>> ma = MiniArray(1)
>>> beam = ma.beam(
        sky=sky,
        pointing=ma_pointing
    )

Calling print() on a Sky object enables the display of its value attribute structure (which matches the definition of the sky instance):

>>> print(beam)
<class 'nenupy.astro.sky.HpxSky'> instance
value: (1, 3, 2, 49152)
    * time: (1,)
    * frequency: (3,)
    * polarization: (2,)
    * coordinates: (49152,)

To plot() the computed Mini-Array response at 75 MHz, in NE polarization:

>>> beam[0, 2, 1].plot(
        decibel=True,
        colorbar_label=''
    )
_images/ma1_beam.png
beamsquint_correction(coords, frequency=<Quantity 50. MHz>)[source]

Corrects for the beamsquint effect.

Example:
>>> from astropy.coordinates import SkyCoord, AltAz
>>> from astropy.time import Time
>>> import astropy.units as u
>>> from nenupy import nenufar_position
>>> from nenupy.instru import MiniArray
>>> position = SkyCoord(
        0*u.deg,
        30*u.deg,
        frame=AltAz(
            obstime=Time("2021-01-01 12:00:00"),
            location=nenufar_position
        )
    )
>>> ma = MiniArray()
>>> corrected_position = ma.beamsquint_correction(
        coords=position,
        frequency=50*u.MHz
    )
>>> corrected_position.az.deg, corrected_position.alt.deg
(0., 22.91422672)
effective_area(frequency=<Quantity 50. MHz>, elevation=<Quantity 90. deg>)[source]

Computes the effective area of a NenuFAR Mini-Array. The effective area of a Mini-Array (\(\mathcal{A}_{\rm eff,\ MA}\)) is computed as the sum of dipole effective areas (\(\mathcal{A}_{\rm eff, ant}\)), while taking into account overlaps. This is a function of frequency (\(\nu\)) and elevation (\(\theta\)):

\[\mathcal{A}_{\rm eff,\ MA} (\nu) = \sum_{\rm ant} \mathcal{A}_{\rm eff, ant} (\nu) \sin( \theta )\]

with

\[\mathcal{A}_{\rm eff, ant} (\nu) = \frac{\lambda^2}{3}\]

the NenuFAR dipole antenna effective area.

Parameters:
  • frequency (Quantity) – Frequency at which the effective area is computed. Default is 50 MHz.

  • elevation (Quantity) – Elevation at which the effective area is computed. Default is 90 deg, i.e., as seen from the zenith.

Returns:

Effective area of a Mini-Array shaped as frequency.

Return type:

Quantity

Example:
>>> from nenupy.instru import MiniArray
>>> import astropy.units as u
>>> ma = MiniArray()
>>> ma.effective_area(50*u.MHz)
227.68377 m2
>>> ma = MiniArray()
>>> ma.effective_area(frequency=50*u.MHz, elevation=45*u.deg)
160.99673 m2
>>> ma = MiniArray()["Ant01"]
>>> ma.effective_area(50*u.MHz)
11.979179 m2
>>> ma = MiniArray()
>>> ma.effective_area(u.Quantity([20, 30, 40], unit='MHz'))
[693.44216, 532.97815, 355.85306] m2

See also

Effective area

property index

Mini-Array index. ‘Core’ Mini-Arrays have indices ranging from 0 to 95. ‘Remote’ Mini-Arrays have indices ranging from 100 to 105.

Setter:

Mini-Array index.

Getter:

Mini-Array index.

Type:

int

static instrument_temperature(frequency=<Quantity 50. MHz>, lna_filter=0)[source]

Instrument temperature at a given frequency. This depends on the Low Noise Amplifier characteristics.

Parameters:
  • frequency (Quantity) – Frequency at which computing the instrument temperature. Default is 50 MHz.

  • lna_filter (int) – Local Noise Amplifier high-pass filter selection. Available values are 0, 1, 2, 3. They correspond to minimal frequencies 10, 15, 20, 25 MHz respectively. Default is 0, i.e., 10 MHz filter.

Returns:

Instrument temperature in Kelvins

Return type:

Quantity

Warning

For the time being, only lna_filter values 0 and 3 are available.

Example:
>>> from nenupy.instru import MiniArray
>>> import astropy.units as u
>>> ma = MiniArray()
>>> ma.instrument_temperature(frequency=70*u.MHz)
526.11213 K
property rotation

Mini-Array rotation. Each NenuFAR Mini-Array has its own rotation with respect to the others by angles multiple of 10 deg.

Setter:

Mini-Array rotation.

Getter:

Mini-Array rotation.

Type:

Quantity

class nenupy.instru.nenufar.NenuFAR(miniarray_antennas=array([0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17, 18]), include_remote_mas=False)[source]

Bases: Interferometer

Main class to handle a NenuFAR array.

Added in version 2.0.0.

Parameters:
  • miniarray_antennas (numpy.ndarray or slice) – Mini-Arrays antennas selection. Default is numpy.r_[:19], i.e., the full 19 dipole antennas. See MiniArray for different input values.

  • include_remote_mas (bool) – Include or not the remote Mini-Arrays. Default is False, i.e., only the dense ‘core’ of 96 Mini-Arrays is considered.

Example:

Instantiating NenuFAR:

>>> from nenupy.instru import NenuFAR
>>> nenufar = NenuFAR()

Sub-arraying on an existing NenuFAR instance:

>>> sub_nenufar = NenuFAR()["MA001", "MA002", "MA104"]
>>> sub_nenufar.antenna_names
array(['MA001', 'MA002'], dtype='<U5')

If include_remote_mas is True, the remote Mini-Arrays are included in the array and selecting MA104 as above would take this remote Mini-Array into account:

>>> sub_nenufar = NenuFAR(include_remote_mas=True)["MA001", "MA002", "MA104"]
>>> sub_nenufar.antenna_names
array(['MA001', 'MA002', 'MA104'], dtype='<U5')

Combining two NenuFAR instances:

>>> nenufar1 = NenuFAR()["MA001", "MA006"]
>>> nenufar2 = NenuFAR()["MA010", "MA056"]
>>> resulting_array = nenufar1 + nenufar2
>>> resulting_array.antenna_names
array(['MA001', 'MA006', 'MA010', 'MA056'], dtype='<U5')

Note

The result of the addition operation, namely resulting_array in this example will conserve the properties of the first member, namely nenufar1. This is particularly true for the attributes include_remote_mas and miniarray_antennas.

See also

More details on this class usage can be found in Array Configuration and Instrument Properties.

Attributes Summary

miniarray_antennas

List of Mini-Array antennas.

include_remote_mas

miniarray_rotations

position

Array's position.

antenna_names

Antenna names.

antenna_positions

Antenna positions.

antenna_gains

Antenna gains.

baselines

Instrument baselines.

size

Number of elements belonging to the array.

Methods Summary

beam(sky, pointing[, analog_pointing, ...])

Computes the NenuFAR beam over the sky for a given pointing.

effective_area([frequency, elevation])

Computes the effective area of NenuFAR.

instrument_temperature([frequency, lna_filter])

Instrument temperature at a given frequency.

plot(**kwargs)

Plots the antenna distribution.

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

Computes the array factor of the antenna distribution.

system_temperature([frequency, ...])

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

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.

angular_resolution([frequency])

Computes the angular resolution of the antenna array.

confusion_noise([frequency, lofar])

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

Attributes and Methods Documentation

beam(sky, pointing, analog_pointing=None, configuration=<nenupy.instru.nenufar.NenuFAR_Configuration object>, return_complex=False, normalize=True)[source]

Computes the NenuFAR beam over the sky for a given pointing.

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

where \(\nu\) is the frequency, \(\phi\) is the azimuth, \(\theta\) is the elevation, \(\mathcal{G}_{\rm MA}\) is the MiniArray response (see beam()) and \(\mathcal{F}_{\rm NenuFAR}\) 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. (time, frequency, polarization and coordinates are used as inputs for the computation).

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

  • analog_pointing (Pointing) – Instance of Pointing that defines the analog pointing directions over the time. This pointing is subject to beamsquint corrections.

  • configuration (NenuFAR_Configuration) – NenuFAR configuration to consider during the beam simulation. The beamsquint correction and its frequency setting are defined here. Default is NenuFAR_Configuration(beamsquint_correction=True, beamsquint_frequency=50MHz).

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

effective_area(frequency=<Quantity 50. MHz>, elevation=<Quantity 90. deg>)[source]

Computes the effective area of NenuFAR. The effective area of NenuFAR (\(\mathcal{A}_{\rm eff,\ NenuFAR}\)) is computed as \(n_{\rm Mini-Arrays}\) times the effective area of one Mini-Array (\(\mathcal{A}_{\rm eff,\ MA}\)) as a function of the frequency \(\nu\), where \(n_{\rm Mini-Arrays}\) is the number of Mini-Arrays included. This method also takes into account the active antennas within each Mini-Array (such as defined by miniarray_antennas).

\[\mathcal{A}_{\rm eff,\ NenuFAR} (\nu) = n_{\rm Mini-Arrays} \mathcal{A}_{\rm eff,\ MA} (\nu)\]
Parameters:
  • frequency (Quantity) – Frequency at which the effective area is computed. Default is 50 MHz.

  • elevation (Quantity) – Elevation at which the effective area is computed. Default is 90 deg, i.e., as seen from the zenith.

Returns:

Effective area of NenuFAR shaped as frequency.

Return type:

Quantity

Example:
>>> from nenupy.instru import NenuFAR
>>> import astropy.units as u
>>> nenufar = NenuFAR()
>>> nenufar.effective_area(50*u.MHz)
18214.701 m2
>>> from nenupy.instru import NenuFAR
>>> import astropy.units as u
>>> nenufar = NenuFAR()
>>> nenufar.effective_area(u.Quantity([20, 30, 40], unit='MHz'))
[55475.372, 42638.252, 28468.245] m2

See also

effective_area() for the computation of \(\mathcal{A}_{\rm eff,\ MA}\) and Effective area.

property include_remote_mas
static instrument_temperature(frequency=<Quantity 50. MHz>, lna_filter=0)[source]

Instrument temperature at a given frequency. This depends on the Low Noise Amplifier characteristics.

Parameters:
  • frequency (Quantity) – Frequency at which computing the instrument temperature. Default is 50 MHz.

  • lna_filter (int) – Local Noise Amplifier high-pass filter selection. Available values are 0, 1, 2, 3. They correspond to minimal frequencies 10, 15, 20, 25 MHz respectively. Default is 0, i.e., 10 MHz filter.

Returns:

Instrument temperature in Kelvins

Return type:

Quantity

Example:
>>> from nenupy.instru import MiniArray
>>> import astropy.units as u
>>> ma = MiniArray()
>>> ma.instrument_temperature(frequency=70*u.MHz)
526.11213 K
property miniarray_antennas

List of Mini-Array antennas.

property miniarray_rotations
class nenupy.instru.nenufar.NenuFAR_Configuration(beamsquint_correction=True, beamsquint_frequency=<Quantity 50. MHz>)[source]

Bases: object

property beamsquint_frequency
class nenupy.instru.nenufar.Polarization(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)[source]

Bases: Enum

Enumerator of the different available polarizations of NenuFAR.

NE = <nenupy.instru.nenufar._AntennaGain object>
NW = <nenupy.instru.nenufar._AntennaGain object>