nenupy.io.xst.TV_Nearfield

class nenupy.io.xst.TV_Nearfield(nearfield, source_imprints, npix, time, frequency, radius, mini_arrays, stokes)[source]

Bases: object

Class to handle NenuFAR TV nearfield storage and display.

__init__(nearfield, source_imprints, npix, time, frequency, radius, mini_arrays, stokes)[source]

Produce an instance of TV_Nearfield.

Parameters:
  • nearfield (ndarray) – The Near-field map to display, image must be square

  • source_imprints (dict) – Dictionnary of Near-field imprints from radio sources, each imprint must be of same shape as the main near-field image

  • npix (int) – Side pixel width of the nearfield image

  • time (Time) – Average time of the near-field data acquisition

  • frequency (Quantity) – Average frequency of the near-field data acquisition

  • radius (Quantity) – Ground radius of projected near-field

  • mini_arrays (ndarray) – Mini-Arrays indices used to compute the near-field

  • stokes (str) – Stokes parameter represented

Raises:

ValueError – Raised if nearfield is not square, or source_imprints do not match nearfield, or npix does not match the nearfield dimensions

Example

>>> from nenupy.io.xst import XST
>>> from nenupy.io.xst import TV_Nearfield
>>> xst = XST(".../nenupy/tests/test_data/XST.fits")
>>> xx_data = xst.get(polarization="XX")
>>> nearfield, src_dict = xx_data.make_nearfield(
        radius=400*u.m,
        npix=64,
        sources=["Tau A"]
    )
>>> tv_nf = TV_Nearfield(
        nearfield=nearfield,
        source_imprints=src_dict,
        npix=nearfield.shape[0],
        time=xst.time[0],
        frequency=xst.frequencies[0],
        radius=400*u.m,
        mini_arrays=xst.mini_arrays,
        stokes="XX",
    )

Methods

__init__(nearfield, source_imprints, npix, ...)

Produce an instance of TV_Nearfield.

from_fits(file_name)

Load a nearfield previously stored in a FITS file.

save_fits(file_name)

Save a near-field made from NenuFAR-TV data as a FITS file.

save_png([figname, fig, ax, figsize, ...])

Display the near-field image.

classmethod from_fits(file_name)[source]

Load a nearfield previously stored in a FITS file.

Parameters:

file_name (str) – Path to the FITS file containing a near-field image (whose format is such as created by the save_fits() method).

Returns:

Instance of TV_Nearfield.

Return type:

TV_Nearfield

Example

>>> from nenupy.io.xst import TV_Nearfield
>>> nf = TV_Nearfield.from_fits("/path/to/nearfield.fits")
save_fits(file_name)[source]

Save a near-field made from NenuFAR-TV data as a FITS file.

Parameters:

file_name (str) – Name of the file to save

Example

>>> from nenupy.io.xst import NenufarTV
>>> tv = NenufarTV("20191204_132113_nenufarTV.dat")
>>> nf_object = tv.compute_nearfield_tv(sources=["Cyg A", "Cas A", "Sun"])
>>> nf_object.save_fits(file_name="/path/to/nearfield.fits")
save_png(figname='', fig=None, ax=None, figsize=(10, 10), decibel=True, cmap='YlGnBu_r', vmin=None, vmax=None, cbar_label=None, title=None)[source]

Display the near-field image.

Parameters:
  • figname (str, optional) – Name of the figure, if given the figure will be saved, by default “”

  • fig (Figure, optional) – matplotlib figure instance, by default None

  • ax (Axes, optional) – matplotlib ax instance, by default None

  • figsize (Tuple[int, int], optional) – Size of the figure, by default (10, 10)

  • decibel (bool, optional) – Set the scale of the data displayed to dB (i.e., \({\rm dB} = 10 \log_{10}({\rm data})\)), by default True

  • cmap (str, optional) – Color map used to represent the data, by default “YlGnBu_r”

  • vmin (float, optional) – Minimal value displayed, by default None (i.e., overall minimal value)

  • vmax (float, optional) – Maximal value displayed, by default None (i.e., overall maximal value)

  • cbar_label (str, optional) – Label of the color bar, by default None (i.e., automatic label)

  • title (str, optional) – Title of the plot, by default None (i.e., automatic label)

Example

>>> from nenupy.io.xst import XST, TV_Nearfield

>>> xst = XST(".../nenupy/tests/test_data/XST.fits")
>>> xx_data = xst.get(polarization="XX")
>>> nearfield, src_dict = xx_data.make_nearfield(
        radius=400*u.m,
        npix=64,
        sources=["Tau A"]
    )
>>> tv_nf = TV_Nearfield(
        nearfield=nearfield,
        source_imprints=src_dict,
        npix=nearfield.shape[0],
        time=xst.time[0],
        frequency=xst.frequencies[0],
        radius=400*u.m,
        mini_arrays=xst.mini_arrays,
        stokes="XX",
    )
>>> tv_nf.save_png(figname="")
../_images/make_nearfield_xst.png