nenupy.instru.instrument_tools


NenuFAR Tools

nenupy.instru.instrument_tools.freq2sb(frequency)[source]

Conversion between the frequency \(\nu\) and the NenuFAR sub-band index \(n_{\rm SB}\). Each NenuFAR sub-band has a bandwidth of \(\Delta \nu = 195.3125\, \rm{kHz}\):

\[n_{\rm SB} = \lfloor*{ \frac{\nu}{\Delta \nu} + \frac{1}{2} \rfloor\]
Parameters:

frequency (Quantity) – Frequency to convert in sub-band index.

Returns:

Sub-band index, same shape as frequency.

Return type:

int or ndarray

Example:
from nenupy.instru import freq2sb
import astropy.units as u

freq2sb(frequency=50.5*u.MHz)
freq2sb(frequency=[50.5, 51]*u.MHz)
nenupy.instru.instrument_tools.generate_nenufar_subarrays(n_subarrays=2, include_remote_mas=False)[source]

Generates NenuFAR sub-arrays of Mini-Arrays. The sub-arraying is done completely randomly.

Parameters:
  • n_subarrays (int) – Number of sub-arrays to generate from the NenuFAR Mini-Array distribution. Default is 2.

  • include_remote_mas (bool) – Include or not the remote Mini-Arrays.

Returns:

A list of Mini-Array names for each sub-array.

Return type:

list of ndarray

Example:
from nenupy.instru import generate_nenufar_subarrays
from nenupy.instru import NenuFAR

sub_arrays = generate_nenufar_subarrays(n_subarrays=2)
sub_array_1 = NenuFAR()[sub_arrays[0]]
sub_array_2 = NenuFAR()[sub_arrays[1]]
nenupy.instru.instrument_tools.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 instrument_temperature
import astropy.units as u

instrument_temperature(frequency=70*u.MHz)
nenupy.instru.instrument_tools.miniarrays_rotated_like(rotations=[0])[source]

Returns the Mini-Array indices whose rotations match the rotations argument. A \(60^{\circ}\) modulo is automatically applied to all rotation parameters.

Parameters:

rotations (list`[`int]) – Mini-Array rotation(s) to select. A ValueError is raised if the values are not integers and/or if they are not multiples of 10.

Returns:

Mini-Array indices.

Return type:

ndarray

Example:
from nenupy.instru import miniarrays_rotated_like

miniarrays_rotated_like([10])
nenupy.instru.instrument_tools.read_cal_table(calibration_file=None)[source]

Reads NenuFAR antenna delays calibration file.

Parameters:

calibration_file (str) – Name of the calibration file to read. If None or 'default' the standard calibration file is read.

Returns:

Antenna delays shaped as (frequency, mini-arrays, polarizations).

Return type:

ndarray

nenupy.instru.instrument_tools.sb2freq(subband)[source]

Conversion between NenuFAR sub-band index \(n_{\rm SB}\) to sub-band starting frequency \(\nu_{\rm start}\).

\[\nu_{\rm start} = n_{\rm SB} \times \Delta \nu\]

Each NenuFAR sub-band has a bandwidth of \(\Delta \nu = 195.3125\, \rm{kHz}\), therefore, the sub-band \(n_{\rm SB}\) goes from \(\nu_{\rm start}\) to \(\nu_{\rm stop} = \nu_{\rm start} + \Delta \nu\).

Parameters:

subband (int or ndarray of int) – Sub-band index (from 0 to 511).

Returns:

Sub-band start frequency \(\nu_{\rm start}\) in MHz.

Return type:

Quantity

Example:
from nenupy.instru import sb2freq

sb2freq(subband=1)
sb2freq(subband=[1, 2, 3, 4])