nenupy.beamlet.bstdata
BST_Data
Object class BST_Data
(from the beamlet
module)
inherits from Beamlet
.
It is designed to read the so called Beamformed Statistics Data
or BST from the
NenuFAR
low-frequency radio-telescope.
Initialization
The BST_Data
object needs to be provided with a path to
a BST observation file (in FITS format):
>>> from nenupy.beamlet import BST_Data
>>> bst = BST_Data(
bstfile='/path/to/BST.fits'
)
Data selection
Data selection is enabled thanks to dedicated keywords
that set BST_Data
attributes. They can be set
directly or passed as keyword arguments to
BST_Data.select()
:
BST_Data.abeam
: Analog beam index selection.BST_Data.dbeam
: Digital beam index selection.BST_Data.polar
: Polarization selection.BST_Data.timerange
: Time range selection.BST_Data.freqrange
: Frequency range selection.
>>> from nenupy.beamlet import BST_Data
>>> bst = BST_Data(
bstfile='/path/to/BST.fits'
)
>>> data = bst.select(
timerange=['2020-01-01 12:00:00', '2020-01-01 13:00:00'],
freqrange=[55, 65],
polar='NW',
dbeam=0
)
Observation and current selection informations
Once initialized, a BST_Data
contains many
informations regarding the observation accessible via specific
getters such as:
BST_Data.abeams
: Available analog beal indicesBST_Data.dbeams
: Available digital beam indicesBST_Data.mas
: Mini-Arrays per perBST_Data.abeam
BST_Data.marot
: Mini-Array rotations perBST_Data.abeam
BST_Data.ants
: Antennas used perBST_Data.abeam
BST_Data.beamlets
: Beamlet indices perBST_Data.dbeam
BST_Data.freqs
: Available frequencies perBST_Data.dbeam
BST_Data.freq
: Frequency selectionBST_Data.time
: Time selectionBST_Data.azana
: Azimuth pointed perBST_Data.abeam
BST_Data.elana
: Elevation pointed perBST_Data.abeam
BST_Data.utana
: Time range of analog pointings.BST_Data.azdig
: Azimuth pointed perBST_Data.dbeam
BST_Data.eldig
: Elevation pointed perBST_Data.dbeam
BST_Data.utdig
: Time range of digital pointings.
Other informations linked to the base class
nenupy.beamlet.beamlet.Beamlet
can also be accessed.
- class nenupy.beamlet.bstdata.BST_Data(bstfile)[source]
Bases:
Beamlet
Class to read NenuFAR BST data stored as FITS files.
- Parameters:
bstfile (str) – Path to BST file.
- property abeam
Analog beam index selected for the available beams (see
BST_Data.abeams
).- Setter:
Analog beam index
- Getter:
Analog beam index
- Type:
- property abeams
Analog beam indices recorded during the observation. Each analog beam defines a particular set of Mini-Arrays as well as a given pointing direction.
- Getter:
Analog beam indices available
- Type:
- property ants
Antenna used within a Mini-Array. This is analog beam dependant (
BST_Data.abeam
).- Getter:
Antenna selected
- Type:
- property azana
Pointed Azimuth for a given analog beam.
- Getter:
Azimuth
- Type:
- property azdig
Pointed Azimuth for a given digital beam.
- Getter:
Azimuth
- Type:
- property beamlets
Beamlet indices associated with a given digital beam (set with
BST_Data.dbeam
).- Getter:
Beamlet indices
- Type:
- property bstfile
Path toward BST FITS file.
- Setter:
BST file
- Getter:
BST file
- Type:
- Example:
>>> from nenupy.beamlet import BST_Data >>> bst = BST_Data( bstfile='/path/to/BST.fits' ) >>> bst.bstfile '/path/to/BST.fits'
- property dbeam
Digital beam index selected for the available beams (see
BST_Data.dbeams
).- Setter:
Digital beam index
- Getter:
Digital beam index
- Type:
- property dbeams
Digital beam indices recorded during the observation. Each digital beam corresponds to a given pointing direction. They are also linked to corresponding
BST_Data.abeams
.- Getter:
Digital beam indices available
- Type:
- property elana
Pointed Elevation for a given analog beam.
- Getter:
Elevation
- Type:
- property eldig
Pointed Elevation for a given digital beam.
- Getter:
Elevation
- Type:
- property freq
Current frequency selection made from setting
BST_Data.freqrange
.- Getter:
Frequency selection
- Type:
- property freqrange
Frequency selection.
- Setter:
Frequency range
- Getter:
Frequency range
- Type:
- Example:
Set frequency to the closest available from 55 MHz:
>>> bst.freqrange = 55 >>> bst.freq <Quantity [54.98047] MHz>
Select frequencies falling between two values: (assumed to be defined in MHz):
>>> bst.freqrange = [55, 56] >>> bst.freq <Quantity [55.17578 , 55.371094, 55.566406, 55.76172 , 55.95703 ] MHz>
or, using
astropy.units.Quantity
:>>> import astropy.units as u >>> import numpy as np >>> bst.freqrange = np.array([55e6, 56e6]) * u.Hz >>> bst.freq <Quantity [55.17578 , 55.371094, 55.566406, 55.76172 , 55.95703 ] MHz>
- property freqs
Recorded frequencies for a given digital beam (set with
BST_Data.dbeam
). Converted to mid sub-band frequencies- Getter:
Available frequencies
- Type:
- property marot
Rotations corresponding to Mini-Arrays listed in
BST_Data.mas
.- Getter:
Mini-Arrays rotations
- Type:
- property mas
Mini-Arrays that have been used for a given analog beam (set with
BST_Data.abeam
).- Getter:
Mini-Arrays selected
- Type:
- property polar
Polarization selection (between
'NW'
and'NE'
).- Setter:
Polarization
- Getter:
Polarization
- Type:
- select(**kwargs)[source]
Select data from the BST observation file.
- Parameters:
freqrange (
float
,list
,numpy.ndarray
orastropy.units.Quantity
, optional) – Frequency range (seeBST_Data.freqrange
).timerange (
str
,list
,numpy.ndarray
orastropy.time.Time
) – Time range (seeBST_Data.timerange
).dbeam (
int
, optional) – Digital beam index (seeBST_Data.dbeam
)polar (
str
, optional) – Polarization (seeBST_Data.polar
)
- Returns:
Selected data
- Return type:
- property time
Current time selection made from setting
BST_Data.timerange
.- Getter:
Time selection
- Type:
- property timerange
Time selection.
- Setter:
Time range
- Getter:
Time range
- Type:
- Example:
Set time to the closest available from a given ISOT value:
>>> bst.timerange = '2019-11-29T15:00:00' >>> bst.time <Time object: scale='utc' format='jd' value=[2458817.125]>
Select times falling between two values:
>>> bst.timerange = ['2019-11-29T15:00:00', '2019-11-29T16:00:00'] >>> bst.time <Time object: scale='utc' format='jd' value=[2458817.125 ... 2458817.16666667]>
or, using
astropy.time.Time
:>>> from astropy.time import Time >>> t0 = Time('2019-11-29T15:00:00') >>> t1 = Time('2019-11-29T16:00:00') >>> bst.timerange = [t0, t1] >>> bst.time <Time object: scale='utc' format='jd' value=[2458817.125 ... 2458817.16666667]>
or:
>>> from astropy.time import Time, TimeDelta >>> t0 = Time('2019-11-29T15:00:00') >>> t1 = t0 + TimeDelta(3600, format='sec') >>> bst.timerange = [t0, t1] >>> bst.time <Time object: scale='utc' format='jd' value=[2458817.125 ... 2458817.16666667]>
- property utana
This returns the time range of all pointings associated with a given
BST_Data.abeam
. The shape of the output array is:[[pointing_1_start, pointing_1_stop], [pointing_2_start, pointing_2_stop], …]
- Getter:
Time range of analog pointings.
- Type:
- property utdig
This returns the time range of all pointings associated with a given
BST_Data.dbeam
. The shape of the output array is:[[pointing_1_start, pointing_1_stop], [pointing_2_start, pointing_2_stop], …]
- Getter:
Time range of digital pointings.
- Type: