nenupy.io.tf.Spectra

class nenupy.io.tf.Spectra(filename, check_missing_data=True)[source]

Bases: object

Class to read UnDySPuTeD Time-Frequency files (.spectra extension).

Notes

Added in version 2.6.0.

__init__(filename, check_missing_data=True)[source]

Generate an instance of Spectra.

Parameters:
  • filename (str) – Path to the UnDySPuTeD file to read.

  • check_missing_data (bool, optional) – Check for missing data, not checking may speed up the reading time, by default False

Raises:

ValueError – If the filename extension differs from ‘.spectra’.

Methods

__init__(filename[, check_missing_data])

Generate an instance of Spectra.

get([file_name])

Perform data selection and pipeline computation.

info()

Display informations about the file.

select_raw_data(tmin_unix, tmax_unix, ...)

Select a subset from a time-frequency NenuFAR dataset.

Attributes

filename

Path to the .spectra file.

frequency_max

Highest recorded frequency.

frequency_min

Lowest recorded frequency.

time_max

Final time of the data content.

time_min

Starting time of the data content.

property filename

Path to the .spectra file.

Type:

str

property frequency_max

Highest recorded frequency. This value is defined at the channel granularity, i.e. it corresponds to the last channel of the highest sub-band.

Example

>>> from nenupy.io.tf import Spectra
>>> sp = Spectra("/my/file.spectra")
>>> sp.frequency_max
57.421875 MHz
Type:

Quantity

property frequency_min

Lowest recorded frequency. This value is defined at the channel granularity, i.e. it corresponds to the first channel of the lowest sub-band.

Example

>>> from nenupy.io.tf import Spectra
>>> sp = Spectra("/my/file.spectra")
>>> sp.frequency_min
19.921875 MHz
Type:

Quantity

get(file_name=None, **pipeline_kwargs)[source]

Perform data selection and pipeline computation.

Parameters:
  • file_name (str, default: None) – If different than None (default value), name of the HDF5 file (extension ‘.hdf5’) to create and store the result.

  • **pipeline_kwargs – Any pipeline parameter passed as keyword argument from the list below. Changes applied here are not kept once the method has resolved.

  • tmin (str or Time, default: \({\rm min}(t)\)) – Lower edge of time selection, can either be given as a Time object or an ISOT/ISO string.

  • tmax (str or Time, default: \({\rm max}(t)\)) – Upper edge of time selection, can either be given as an Time object or an ISOT/ISO string.

  • fmin (float or Quantity, default: \({\rm min}(\nu)\)) – Lower frequency boundary selection, can either be given as a Quantity object or float (assumed to be in MHz in that case).

  • fmax (float or Quantity, default: \({\rm max}(\nu)\)) – Higher frequency boundary selection, can either be given as a Quantity object or float (assumed to be in MHz in that case).

  • beam (int, default: first recorded beam) – Beam selection, a single integer corresponding to the index of a recorded numerical beam is expected.

  • dispersion_measure (float or Quantity, default: None) – Enable de-dispersion of the data by this Dispersion Measure. Note that the de_disperse() task should be present in the planned pipeline (pipeline). It can either be provided as a Quantity object or a float (assumed to be in \({\rm pc}\,{\rm cm}^{-3}\) in that case).

  • rotation_measure (float or Quantity, default: None) – Enable the correction of the Faraday rotation using this Rotation Measure. Note that the correct_faraday_rotation() task should be present in the planned pipeline (pipeline). It can either be provided as a Quantity object or a float (assumed to be in \({\rm rad}\,{\rm m}^{-2}\) in that case).

  • rebin_dt (float or Quantity, default: None) – Desired rebinning time resolution, can either be given as a Quantity object or a float (assumed to be in sec in that case). Note that the time_rebin() task should be present in the planned pipeline (pipeline).

  • rebin_df (float or Quantity, default: None) – Desired rebinning frequency resolution, can either be given as a Quantity object or float (assumed to be in kHz in that case). Note that the frequency_rebin() task should be present in the planned pipeline (pipeline).

  • smooth_frequency_profile (bool, default False) – Keyword used if flatten_subband() is in the pipeline: smooth the frequency profile, this option mau result in unexpected behavior if the subbands are not contiguous or the selected beamlets do not belong to the same digital beam.

  • remove_channels (list or ndarray, default: None) – List of subband channels to remove, e.g. remove_channels=[0,1,-1] would remove the first, second (low-freq) and last channels from each subband. Note that the remove_channels() task should be present in the planned pipeline (pipeline).

  • skycoord (SkyCoord, default: None) – Tracked celestial coordinates used for beam and polarization corrections, a SkyCoord object is expected. Note that the correct_polarization() or correct_parallactic_rotation() task should be present in the planned pipeline (pipeline).

  • calib_dt (float or Quantity, default: None) – Time resolution used for beam and polarization corrections, a Quantity object or a float (assumed in seconds) are expected. Note that the correct_polarization() or correct_parallactic_rotation() task should be present in the planned pipeline (pipeline).

  • dreambeam_parallactic (bool, default: True) – DreamBeam parallactic angle correction (along with 'skycoord' and 'calib_dt'), a boolean is expected. Note that the correct_polarization() task should be present in the planned pipeline (pipeline).

  • stokes (str or list[str], default: "I") – Stokes parameter selection, can either be given as a string or a list of strings, e.g. ['I', 'Q', 'V/I']. Note that the get_stokes() task should be present in the planned pipeline (pipeline).

  • ignore_volume_warning (bool, default: False) – Ignore or not (default value) the limit regarding output data volume.

Returns:

Processed data selection.

Return type:

SData

Examples

>>> from nenupy.io.tf import Spectra
>>> sp = Spectra("/my/file.spectra")
>>> result = sp.get(tmin="2024-01-01 12:00:00", tmax="2024-01-01 13:00:00")
info()[source]

Display informations about the file.

Example

>>> from nenupy.io.tf import Spectra

>>> sp = Spectra("/my/file.spectra")
>>> sp.info()
filename: /my/file.spectra
time_min: 2023-05-27T08:39:02.0000050
time_max: 2023-05-27T08:59:34.2445748
dt: 20.97152 ms
frequency_min: 19.921875 MHz
frequency_max: 57.421875 MHz
df: 3.0517578125 kHz
Available beam indices: ['0']
select_raw_data(tmin_unix, tmax_unix, fmin_hz, fmax_hz, beam)[source]

Select a subset from a time-frequency NenuFAR dataset.

Parameters:
  • tmin_unix (float) – Start time of the selection in UNIX format.

  • tmax_unix (float) – Stop time of the selection in UNIX format.

  • fmin_hz (float) – Minimal frequency of the selection in Hz.

  • fmax_hz (float) – Maximal frequency of the selection in Hz.

  • beam (int) – Beam index of the selection.

Raises:

KeyError – If beam does not correspond to any recorded value.

Returns:

Length-3 tuple, respectively containing the frequency in Hz and the time in unix as numpy arrays, then the selected dataset. The latter should be shaped as (time, frequency, 2, 2) where the last two dimensions are the Jones matrix of the electric field cross correlations. If the chosen time arguments lead to an empty selection, a length-3 tuple of None is returned. If the frequency selection is off, the closest subband is returned.

Return type:

(ndarray, ndarray, Array)

Notes

Added in version 2.6.0.

Example

>>> import nenupy
>>> sp = Spectra("/my/file.spectra")
>>> data = sp.select_raw_data(
        tmin_unix=1685176742.000005,
        tmax_unix=1685176802.000005,
        fmin_hz=21e6,
        fmax_hz=22e6,
        beam=0
    )
>>> data[0].shape, data[1].shape, data[2].shape
((384,), (2856,), (2856, 384, 2, 2))
property time_max

Final time of the data content.

Example

>>> from nenupy.io.tf import Spectra
>>> sp = Spectra("/my/file.spectra")
>>> sp.time_max.isot
'2023-05-27T08:59:34.2445748'
Type:

Time

property time_min

Starting time of the data content.

Example

>>> from nenupy.io.tf import Spectra
>>> sp = Spectra("/my/file.spectra")
>>> sp.time_min.isot
'2023-05-27T08:39:02.0000050'
Type:

Time