nenupy.instru.nenufar.NenuFAR
- 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:
InterferometerMain class to handle a NenuFAR array.
Added in version 2.0.0.
- Parameters:
miniarray_antennas (
numpy.ndarrayorslice) – Mini-Arrays antennas selection. Default isnumpy.r_[:19], i.e., the full 19 dipole antennas. SeeMiniArrayfor different input values.include_remote_mas (
bool) – Include or not the remote Mini-Arrays. Default isFalse, 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
NenuFARinstance:>>> sub_nenufar = NenuFAR()["MA001", "MA002", "MA104"] >>> sub_nenufar.antenna_names array(['MA001', 'MA002'], dtype='<U5')
If
include_remote_masisTrue, the remote Mini-Arrays are included in the array and selectingMA104as 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
NenuFARinstances:>>> 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_arrayin this example will conserve the properties of the first member, namelynenufar1. This is particularly true for the attributesinclude_remote_masandminiarray_antennas.
See also
More details on this class usage can be found in Array Configuration and Instrument Properties.
Attributes Summary
List of Mini-Array antennas.
include_remote_masArray's position.
Antenna names.
Antenna positions.
Antenna gains.
Instrument baselines.
Number of elements belonging to the array.
Methods Summary
beam(sky, pointing[, analog_pointing, ...])Computes the NenuFAR beam over the
skyfor a givenpointing.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
- __init__(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]
Methods
__init__([miniarray_antennas, ...])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[, analog_pointing, ...])Computes the NenuFAR beam over the
skyfor a givenpointing.confusion_noise([frequency, lofar])Confusion rms noise \(\sigma_{\rm c}\) (parameter used for specifying the width of the confusion distribution) computed as:
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.
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
Add delay errors for each antennae They could be cable connection errors during construction, cables of wrong length, ...
Antenna gains.
Antenna names.
Antenna positions.
antenna_weightsInstrument baselines.
include_remote_masList of Mini-Array antennas.
Array's position.
Number of elements belonging to the array.
- angular_resolution(frequency=<Quantity 50. MHz>)
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.
- 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.
- property antenna_names
Antenna names.
- Setter:
Array of antenna names.
- Getter:
Array of antenna names.
- Type:
- 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:
- array_factor(sky, pointing, return_complex=False, normalize=True)
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
skyas the desired output (in terms of time, frequency and sky positions). It evaluates the effective pointing directions for every time step defined inskyregarding thepointinginput.- Parameters:
sky (
Sky) – Desired output contained in aSkyinstance. (time,frequencyandcoordinatesare used as inputs for the computation).pointing (
Pointing) – Instance ofPointingthat defines the targeted pointing directions over the time.return_complex (
bool) – Return complex array factor ifTrueor power ifFalsenormalize (
bool) – Return the normalized array factor. Default isTrue.
- Returns:
Array factor of the antenna distribution shaped as
(time, frequency, 1, coordinates).- Return type:
- 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
skyfor a givenpointing.\[\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
skyas the desired output (in terms of time, frequency, polarization and sky positions). It evaluates the effective pointing directions for every time step defined inskyregarding thepointinginput.- Parameters:
sky (
Sky) – Desired output contained in aSkyinstance. (time,frequency,polarizationandcoordinatesare used as inputs for the computation).pointing (
Pointing) – Instance ofPointingthat defines the targeted numerical pointing directions over the time.analog_pointing (
Pointing) – Instance ofPointingthat 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 isNenuFAR_Configuration(beamsquint_correction=True, beamsquint_frequency=50MHz).
- Returns:
The instance of
Skygiven as input is returned, its attributevalueis updated with the result of the beam computation (stored as anArray) and shaped as(time, frequency, polarization, coordinates).- Return type:
See also
- confusion_noise(frequency=<Quantity 50. MHz>, lofar=True)
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 (
floatorQuantity) – Frequency at which computing the confusion noise. In MHz if no unit is provided. Default is50 MHz.miniarrays (
int,listorndarray) – Mini-Array indices to take into account. Default isNone(all available MAs).lofar – If set to
True(recommended), the confusion noise is estimated using Eq. 6 of van Haarlem et al. (2013).
- Type:
- Returns:
Confusion rms noise in Jy/beam
- Return type:
- Example:
>>> import astropy.units as u >>> <instrument>.confusion_noise( frequency=50*u.MHz )
- 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 byminiarray_antennas).\[\mathcal{A}_{\rm eff,\ NenuFAR} (\nu) = n_{\rm Mini-Arrays} \mathcal{A}_{\rm eff,\ MA} (\nu)\]- Parameters:
- Returns:
Effective area of NenuFAR shaped as
frequency.- Return type:
- 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.
- 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 is50 MHz.lna_filter (
int) – Local Noise Amplifier high-pass filter selection. Available values are0, 1, 2, 3. They correspond to minimal frequencies10, 15, 20, 25 MHzrespectively. Default is0, i.e., 10 MHz filter.
- Returns:
Instrument temperature in Kelvins
- Return type:
- Example:
>>> from nenupy.instru import MiniArray >>> import astropy.units as u >>> ma = MiniArray() >>> ma.instrument_temperature(frequency=70*u.MHz) 526.11213 K
See also
- property miniarray_antennas
List of Mini-Array antennas.
- property miniarray_rotations
- plot(**kwargs)
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 isTrue.patches (
tuple`(`list, colors) of Polygons) – Matplotlib Polygons
- property position
Array’s position.
- Setter:
Position of the array.
- Getter:
Position of the array.
- Type:
- sefd(frequency=<Quantity 50. MHz>, elevation=<Quantity 90. deg>, efficiency=1.0, decoherence=1.0, source_spectrum={}, **kwargs)
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 if50 MHz.elevation (
Quantity) – Pointing elevation impacting theeffective_area(). Default is90 deg.efficiency (
float) – Effective area reducing factor. Default is1., it cannot be greater than1..decoherence (
float) – Parameter that reflects other uncertainties (particularly the unperfect phasing system). Default is1..source_spectrum (
dictofcallable) – 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 acallableobject that takes as inputs a frequency array (of typeQuantity) and returns the source flux density in Jansky (of type (asQuantity).
- Returns:
SEFD in Janskys.
- Return type:
- 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)
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_spectrumargument, seesefd()).- 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 if50 MHz.mode (
ObservingMode) – Observing mode, eitherObservingMode.BEAMFORMINGorObservingMode.IMAGING, default is the former.dt (
Quantity) – Integration time. Default is1 sec.df (
Quantity) – Observing bandwidth. Default is195.3125 kHz.elevation (
Quantity) – Pointing elevation impacting theeffective_area(). Default is90 deg.efficiency (
float) – Effective area reducing factor. Default is1., it cannot be greater than1..decoherence (
float) – Parameter that reflects other uncertainties (particularly the unperfect phasing system). Default is1..source_spectrum (
dictofcallable) – 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 acallableobject that takes as inputs a frequency array (of typeQuantity) and returns the source flux density in Jansky (of typeQuantity).
- Returns:
Array sensitivity.
- Return type:
- 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 )
See also
- system_temperature(frequency=<Quantity 50. MHz>, source_spectrum={}, efficiency=1.0, elevation=<Quantity 90. deg>, **kwargs)
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_spectrumargument 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
efficiencyof the effective area \(A_{\rm eff}\).- Parameters:
frequency (
Quantity) – Frequency for the System Temperature computation. Default is50 MHz.elevation (
Quantity) – Pointing elevation impacting theeffective_area(). Default is90 deg.efficiency (
float) – Effective area reducing factor. Default is1., it cannot be greater than1..source_spectrum (
dictofcallable) – 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 acallableobject that takes as inputs a frequency array (of typeQuantity) and returns the source flux density in Jansky (of typeQuantity).
- Returns:
System Temperature in Kelvins.
- Return type:
See also