Analytic Beams

pyuvdata defines several analytic primary beams for radio telescopes. While these are not realistic models for true antennas (like those represented in pyuvdata.UVBeam), they can be useful in simulation because they are lightweight and fast to evaluate (as opposed to having to interpolate).

The analytic beams defined in pyuvdata are based on a base class, pyuvdata.analytic_beam.AnalyticBeam, which ensures a standard interface and can be used to define other analytic beams in a consistent way (see the analytic beam tutorial). To evaluate analytic beams in particular directions at particular frequencies, use the pyuvdata.analytic_beam.AnalyticBeam.efield_eval() or pyuvdata.analytic_beam.AnalyticBeam.power_eval() methods as appropriate.

The AnalyticBeam base class also provides a yaml constructor that can enable analytic beams to be instantiated directly from yaml files (see yaml constructors, similar constructors are also available for UVBeam objects) and a plugin infrastructure that can automatically include any imported subclass even if they are defined in other packages. This infrastructure, along with the pyuvdata.BeamInterface class, can simplify the setup for simulations.

class pyuvdata.analytic_beam.AnalyticBeam(*, feed_array: ndarray[tuple[Any, ...], dtype[str_]] | list[str] | None = None, feed_angle: ndarray[tuple[Any, ...], dtype[floating]] | list[float] | None = None, mount_type: Literal['alt-az', 'equatorial', 'orbiting', 'x-y', 'alt-az+nasmyth-r', 'alt-az+nasmyth-l', 'phased', 'fixed', 'other'] | None = 'fixed', include_cross_pols: dataclasses.InitVar[bool] = True, x_orientation: dataclasses.InitVar[Optional[Literal['east', 'north']]] = 'east')[source]

Analytic beam base class.

feed_array

Feeds to define this beam for, e.g. x & y or r & l.

Type:

array-like of str

feed_angle

Position angle of a given feed, units of radians. A feed angle of 0 is typically oriented toward zenith for steerable antennas, otherwise toward north for fixed antennas (e.g., HERA, LWA). More details on this can be found on the “Conventions” page of the docs.

Type:

array-like of float

mount_type

Antenna mount type, which describes the optics of the antenna in question. Supported options include: “alt-az” (primary rotates in azimuth and elevation), “equatorial” (primary rotates in hour angle and declination) “orbiting” (antenna is in motion, and its orientation depends on orbital parameters), “x-y” (primary rotates first in the plane connecting east, west, and zenith, and then perpendicular to that plane), “alt-az+nasmyth-r” (“alt-az” mount with a right-handed 90-degree tertiary mirror), “alt-az+nasmyth-l” (“alt-az” mount with a left-handed 90-degree tertiary mirror), “phased” (antenna is “electronically steered” by summing the voltages of multiple elements, e.g. MWA), “fixed” (antenna beam pattern is fixed in azimuth and elevation, e.g., HERA), and “other” (also referred to in some formats as “bizarre”). See the “Conventions” page of the documentation for further details.

Type:

str

Parameters:
  • feed_array (array-like of str) – Feeds to define this beam for, e.g. x & y or r & l.

  • feed_angle (array-like of float) – Position angle of a given feed, units of radians. A feed angle of 0 is typically oriented toward zenith for steerable antennas, otherwise toward north for fixed antennas (e.g., HERA, LWA). More details on this can be found on the “Conventions” page of the docs.

  • mount_type (str) – Antenna mount type, which describes the optics of the antenna in question. Supported options include: “alt-az” (primary rotates in azimuth and elevation), “equatorial” (primary rotates in hour angle and declination) “orbiting” (antenna is in motion, and its orientation depends on orbital parameters), “x-y” (primary rotates first in the plane connecting east, west, and zenith, and then perpendicular to that plane), “alt-az+nasmyth-r” (“alt-az” mount with a right-handed 90-degree tertiary mirror), “alt-az+nasmyth-l” (“alt-az” mount with a left-handed 90-degree tertiary mirror), “phased” (antenna is “electronically steered” by summing the voltages of multiple elements, e.g. MWA), “fixed” (antenna beam pattern is fixed in azimuth and elevation, e.g., HERA), and “other” (also referred to in some formats as “bizarre”). See the “Conventions” page of the documentation for further details. Default is “fixed”.

  • include_cross_pols (bool) – Option to include the cross polarized beams (e.g. xy and yx or en and ne) for the power beam.

  • x_orientation (str) – Physical orientation of the feed for the x feed. Not meaningful for unpolarized analytic beams, but clarifies the orientation of the dipole for for polarized beams like the ShortDipoleBeam and matches with the meaning on UVBeam objects. Used to set values of the feed_angle attribute, ignored in feed_angle is provided. Default is “east”.

get_x_orientation_from_feeds() Literal['east', 'north', None][source]

Get x-orientation equivalent value based on feed information.

Returns:

x_orientation (str) – One of “east”, “north”, or None, based on values present in AnalyticBeam.feed_array and AnalyticBeam.feed_angle.

set_feeds_from_x_orientation(x_orientation)[source]

Set feed information based on x-orientation value.

Populates newer parameters describing feed-orientation (AnalyticBeam.feed_array and AnalyticBeam.feed_angle) based on the “older” x-orientation string. Note that this method will overwrite any previously populated values.

Parameters:

x_orientation (str) – String describing how the x-orientation is oriented. Must be either “north”/ “n”/”ns” (x-polarization of antenna has a position angle of 0 degrees with respect to zenith/north) or “east”/”e”/”ew” (x-polarization of antenna has a position angle of 90 degrees with respect to zenith/north).

validate()[source]

Validate inputs, placeholder for subclasses.

property Naxes_vec

The number of vector axes.

property Nfeeds

The number of feeds.

property Npols

The number of polarizations.

clone(**kw)[source]

Create a new instance of the object with updated parameters.

property east_ind

The index of the east feed in the feed array.

property north_ind

The index of the north feed in the feed array.

efield_eval(*, az_array: ndarray[tuple[Any, ...], dtype[floating]], za_array: ndarray[tuple[Any, ...], dtype[floating]], freq_array: ndarray[tuple[Any, ...], dtype[floating]]) ndarray[tuple[Any, ...], dtype[floating]][source]

Evaluate the efield at the given coordinates.

Parameters:
  • az_array (array_like of floats, optional) – Azimuth values to evaluate the beam at in radians. Must be the same shape as za_array.

  • za_array (array_like of floats, optional) – Zenith values to evaluate the beam at in radians. Must be the same shape as az_array.

  • freq_array (array_like of floats, optional) – Frequency values to evaluate the beam at in Hertz.

Returns:

array_like of complex – An array of beam values. The shape of the evaluated data will be: (Naxes_vec, Nfeeds, freq_array.size, az_array.size)

power_eval(*, az_array: ndarray[tuple[Any, ...], dtype[floating]], za_array: ndarray[tuple[Any, ...], dtype[floating]], freq_array: ndarray[tuple[Any, ...], dtype[floating]]) ndarray[tuple[Any, ...], dtype[floating]][source]

Evaluate the power at the given coordinates.

Parameters:
  • az_array (array_like of floats, optional) – Azimuth values to evaluate the beam at in radians. Must be the same shape as za_array.

  • za_array (array_like of floats, optional) – Zenith values to evaluate the beam at in radians. Must be the same shape as az_array.

  • freq_array (array_like of floats, optional) – Frequency values to evaluate the beam at in Hertz.

Returns:

array_like of float or complex – An array of beam values. The shape of the evaluated data will be: (1, Npols, freq_array.size, az_array.size). The dtype will be a complex type if cross-pols are included, otherwise it will be a float type.

to_uvbeam(freq_array: ndarray[tuple[Any, ...], dtype[floating]], beam_type: Literal['efield', 'power'] = 'efield', pixel_coordinate_system: Literal['az_za', 'orthoslant_zenith', 'healpix'] | None = None, axis1_array: ndarray[tuple[Any, ...], dtype[floating]] | None = None, axis2_array: ndarray[tuple[Any, ...], dtype[floating]] | None = None, nside: int | None = None, healpix_pixel_array: ndarray[tuple[Any, ...], dtype[integer]] | None = None, ordering: Literal['ring', 'nested'] | None = None)[source]

Generate a UVBeam object from an AnalyticBeam object.

This method evaluates the analytic beam at a set of locations and frequencies to create a UVBeam object. This can be useful for testing and some other operations, but it is of course an approximation.

Parameters:
  • freq_array (ndarray of float) – Array of frequencies in Hz to evaluate the beam at.

  • beam_type (str) – Beam type, either “efield” or “power”.

  • pixel_coordinate_system (str) – Pixel coordinate system, options are “az_za”, “orthoslant_zenith” and “healpix”. Forced to be “healpix” if nside is given and by default set to “az_za” if not. Currently, only “az_za” and “healpix” are implemented.

  • axis1_array (ndarray of float) – Coordinates along first pixel axis (e.g. azimuth for an azimuth/zenith angle coordinate system) to evaluate the beam at. Must be regularly spaced. Should not provided for healpix coordinates.

  • axis2_array (ndarray of float) – Coordinates along second pixel axis (e.g. zenith angle for an azimuth/zenith angle coordinate system) to evaluate the beam at. Must be regularly spaced. Should not provided for healpix coordinates.

  • nside (int) – Healpix nside parameter, should only be provided for healpix coordinates.

  • healpix_pixel_array (ndarray of int) – Healpix pixels to include. If nside is provided, defaults to all the pixels in the Healpix map.

  • ordering (str) – Healpix ordering parameter, defaults to “ring” if nside is provided.

Returns:

UVBeam – A new UVBeam object.

plot(*, beam_type: str, freq: float, complex_type: str = 'real', colormap: str | None = None, logcolor: bool | None = None, plt_kwargs: dict | None = None, norm_kwargs: dict | None = None, max_zenith_deg: float = 90.0, savefile: str | None = None)[source]

Make a pretty plot of the beams.

This evaluates the beam on a grid in azimuth and zenith angle and plots the result.

Parameters:
  • freq (int) – The frequency index to plot.

  • complex_type (str) – What to plot for complex beams, options are: [real, imag, abs, phase]. Defaults to “real” for complex beams. Ignored for real beams (i.e. power beams, same feed).

  • colormap (str, optional) – Matplotlib colormap to use. Defaults to “twlight” if complex_type=”phase” and logcolor=False, otherwise it defaults to “viridis” if the data to be plotted is positive definite (e.g. if complex_type=”abs”) and “PRGn” otherwise.

  • logcolor (bool, optional) – Option to use log scaling for the color. Defaults to True for power beams and False for E-field beams. Results in using matplotlib.colors.LogNorm or matplotlib.colors.SymLogNorm if the data have negative values.

  • plt_kwargs (dict, optional) – Keywords to be passed into the matplotlib.pyplot.imshow call.

  • norm_kwargs (dict, optional) – Keywords to be passed into the norm object, typically vmin/vmax, plus linthresh for SymLogNorm.

  • max_zenith_deg (float) – Maximum zenith angle to include in the plot in degrees. Default is 90 to go down to the horizon.

  • savefile (str) – File to save the plot to.

class pyuvdata.AiryBeam(*, feed_array: ndarray[tuple[Any, ...], dtype[str_]] | list[str] | None = None, feed_angle: ndarray[tuple[Any, ...], dtype[floating]] | list[float] | None = None, mount_type: Literal['alt-az', 'equatorial', 'orbiting', 'x-y', 'alt-az+nasmyth-r', 'alt-az+nasmyth-l', 'phased', 'fixed', 'other'] = 'fixed', include_cross_pols: dataclasses.InitVar[bool] = True, x_orientation: dataclasses.InitVar[Optional[Literal['east', 'north']]] = 'east', diameter: float)[source]

A zenith pointed Airy beam.

Airy beams are the diffraction pattern of a circular aperture, so represent an idealized dish. Requires a dish diameter in meters and is inherently chromatic and unpolarized.

The unpolarized nature leads to some results that may be surprising to radio astronomers: if two feeds are specified they will have identical responses and the cross power beam between the two feeds will be identical to the power beam for a single feed.

diameter

Dish diameter in meters.

Type:

float

Parameters:
  • diameter (float) – Dish diameter in meters.

  • include_cross_pols (bool) – Option to include the cross polarized beams (e.g. xy and yx or en and ne) for the power beam.

class pyuvdata.GaussianBeam(*, feed_array: ndarray[tuple[Any, ...], dtype[str_]] | list[str] | None = None, feed_angle: ndarray[tuple[Any, ...], dtype[floating]] | list[float] | None = None, mount_type: Literal['alt-az', 'equatorial', 'orbiting', 'x-y', 'alt-az+nasmyth-r', 'alt-az+nasmyth-l', 'phased', 'fixed', 'other'] = 'fixed', include_cross_pols: dataclasses.InitVar[bool] = True, x_orientation: dataclasses.InitVar[Optional[Literal['east', 'north']]] = 'east', sigma: float | None = None, sigma_type: Literal['efield', 'power'] = 'efield', diameter: float | None = None, spectral_index: float = 0.0, reference_frequency: float = None)[source]

A circular, zenith pointed Gaussian beam.

Requires either a dish diameter in meters or a standard deviation sigma in radians. Gaussian beams specified by a diameter will have their width matched to an Airy beam at each simulated frequency, so are inherently chromatic. For Gaussian beams specified with sigma, the sigma_type defines whether the width specified by sigma specifies the width of the E-Field beam (default) or power beam in zenith angle. If only sigma is specified, the beam is achromatic, optionally both the spectral_index and reference_frequency parameters can be set to generate a chromatic beam with standard deviation defined by a power law:

stddev(f) = sigma * (f/ref_freq)**(spectral_index)

The unpolarized nature leads to some results that may be surprising to radio astronomers: if two feeds are specified they will have identical responses and the cross power beam between the two feeds will be identical to the power beam for a single feed.

sigma

Standard deviation in radians for the gaussian beam. Only one of sigma and diameter should be set.

Type:

float

sigma_type

Either “efield” or “power” to indicate whether the sigma specifies the size of the efield or power beam. Ignored if sigma is None.

Type:

str

diameter

Dish diameter in meters to use to define the size of the gaussian beam, by matching the FWHM of the gaussian to the FWHM of an Airy disk. This will result in a frequency dependent beam. Only one of sigma and diameter should be set.

Type:

float

spectral_index

Option to scale the gaussian beam width as a power law with frequency. If set to anything other than zero, the beam will be frequency dependent and the reference_frequency must be set. Ignored if sigma is None.

Type:

float

reference_frequency

The reference frequency for the beam width power law, required if sigma is not None and spectral_index is not zero. Ignored if sigma is None.

Type:

float

Parameters:
  • sigma (float) – Standard deviation in radians for the gaussian beam. Only one of sigma and diameter should be set.

  • sigma_type (str) – Either “efield” or “power” to indicate whether the sigma specifies the size of the efield or power beam. Ignored if sigma is None.

  • diameter (float) – Dish diameter in meters to use to define the size of the gaussian beam, by matching the FWHM of the gaussian to the FWHM of an Airy disk. This will result in a frequency dependent beam. Only one of sigma and diameter should be set.

  • spectral_index (float) – Option to scale the gaussian beam width as a power law with frequency. If set to anything other than zero, the beam will be frequency dependent and the reference_frequency must be set. Ignored if sigma is None.

  • reference_frequency (float) – The reference frequency for the beam width power law, required if sigma is not None and spectral_index is not zero. Ignored if sigma is None.

  • include_cross_pols (bool) – Option to include the cross polarized beams (e.g. xy and yx or en and ne) for the power beam.

validate()[source]

Post-initialization validation and conversions.

get_sigmas(freq_array: ndarray[tuple[Any, ...], dtype[floating]]) ndarray[tuple[Any, ...], dtype[floating]][source]

Get the sigmas for the gaussian beam using the diameter (if defined).

Parameters:

freq_array (array of floats) – Frequency values to get the sigmas for in Hertz.

Returns:

sigmas (array_like of float) – Beam sigma values as a function of frequency. Size will match the freq_array size.

class pyuvdata.ShortDipoleBeam(*, feed_array: ndarray[tuple[Any, ...], dtype[str_]] | list[str] | None = None, feed_angle: ndarray[tuple[Any, ...], dtype[floating]] | list[float] | None = None, mount_type: Literal['alt-az', 'equatorial', 'orbiting', 'x-y', 'alt-az+nasmyth-r', 'alt-az+nasmyth-l', 'phased', 'fixed', 'other'] | None = 'fixed', include_cross_pols: dataclasses.InitVar[bool] = True, x_orientation: dataclasses.InitVar[Optional[Literal['east', 'north']]] = 'east')[source]

A zenith pointed analytic short dipole beam with two crossed feeds.

A classical short (Hertzian) dipole beam with two crossed feeds aligned east and north. Short dipole beams are intrinsically polarized but achromatic. Does not require any parameters, but the orientation of the dipole labelled as “x” can be specified to align “north” or “east” via the x_orientation parameter (matching the parameter of the same name on UVBeam and UVData objects).

feed_array

Feeds to define this beam for, e.g. x & y.

Type:

array-like of str

feed_angle

Position angle of a given feed, units of radians. A feed angle of 0 is typically oriented toward zenith for steerable antennas, otherwise toward north for fixed antennas (e.g., HERA, LWA). More details on this can be found on the “Conventions” page of the docs.

Type:

array-like of float

mount_type

Antenna mount type, which describes the optics of the antenna in question. Supported options include: “alt-az” (primary rotates in azimuth and elevation), “equatorial” (primary rotates in hour angle and declination) “orbiting” (antenna is in motion, and its orientation depends on orbital parameters), “x-y” (primary rotates first in the plane connecting east, west, and zenith, and then perpendicular to that plane), “alt-az+nasmyth-r” (“alt-az” mount with a right-handed 90-degree tertiary mirror), “alt-az+nasmyth-l” (“alt-az” mount with a left-handed 90-degree tertiary mirror), “phased” (antenna is “electronically steered” by summing the voltages of multiple elements, e.g. MWA), “fixed” (antenna beam pattern is fixed in azimuth and elevation, e.g., HERA), and “other” (also referred to in some formats as “bizarre”). See the “Conventions” page of the documentation for further details.

Type:

str

Parameters:
  • feed_array (array-like of str) – Feeds to define this beam for, e.g. x & y or r & l. Default is [“x”, “y”]

  • feed_angle (array-like of float) – Position angle of a given feed, units of radians. A feed angle of 0 is typically oriented toward zenith for steerable antennas, otherwise toward north for fixed antennas (e.g., HERA, LWA). More details on this can be found on the “Conventions” page of the docs.

  • mount_type (str) – Antenna mount type, which describes the optics of the antenna in question. Supported options include: “alt-az” (primary rotates in azimuth and elevation), “equatorial” (primary rotates in hour angle and declination) “orbiting” (antenna is in motion, and its orientation depends on orbital parameters), “x-y” (primary rotates first in the plane connecting east, west, and zenith, and then perpendicular to that plane), “alt-az+nasmyth-r” (“alt-az” mount with a right-handed 90-degree tertiary mirror), “alt-az+nasmyth-l” (“alt-az” mount with a left-handed 90-degree tertiary mirror), “phased” (antenna is “electronically steered” by summing the voltages of multiple elements, e.g. MWA), “fixed” (antenna beam pattern is fixed in azimuth and elevation, e.g., HERA), and “other” (also referred to in some formats as “bizarre”). See the “Conventions” page of the documentation for further details.

  • include_cross_pols (bool) – Option to include the cross polarized beams (e.g. xy and yx or en and ne) for the power beam.

  • x_orientation (str) – Physical orientation of the feed for the x feed. Not meaningful for unpolarized analytic beams, but clarifies the orientation of the dipole for for polarized beams like the ShortDipoleBeam and matches with the meaning on UVBeam objects. Used to set values of the feed_angle attribute, ignored in feed_angle is provided. Default is “east”.

validate()[source]

Post-initialization validation and conversions.

class pyuvdata.UniformBeam(*, feed_array: ndarray[tuple[Any, ...], dtype[str_]] | list[str] | None = None, feed_angle: ndarray[tuple[Any, ...], dtype[floating]] | list[float] | None = None, mount_type: Literal['alt-az', 'equatorial', 'orbiting', 'x-y', 'alt-az+nasmyth-r', 'alt-az+nasmyth-l', 'phased', 'fixed', 'other'] = 'fixed', include_cross_pols: dataclasses.InitVar[bool] = True, x_orientation: dataclasses.InitVar[Optional[Literal['east', 'north']]] = 'east')[source]

A uniform beam.

Uniform beams have identical responses in all directions, so are quite unphysical but can be useful for testing other aspects of simulators. They are unpolarized and achromatic and do not take any parameters.

The unpolarized nature leads to some results that may be surprising to radio astronomers: if two feeds are specified they will have identical responses and the cross power beam between the two feeds will be identical to the power beam for a single feed.

feed_array

Feeds to define this beam for, e.g. x & y or r & l.

Type:

array-like of str

feed_angle

Position angle of a given feed, units of radians. A feed angle of 0 is typically oriented toward zenith for steerable antennas, otherwise toward north for fixed antennas (e.g., HERA, LWA). More details on this can be found on the “Conventions” page of the docs.

Type:

array-like of float

mount_type

Antenna mount type, which describes the optics of the antenna in question. Supported options include: “alt-az” (primary rotates in azimuth and elevation), “equatorial” (primary rotates in hour angle and declination) “orbiting” (antenna is in motion, and its orientation depends on orbital parameters), “x-y” (primary rotates first in the plane connecting east, west, and zenith, and then perpendicular to that plane), “alt-az+nasmyth-r” (“alt-az” mount with a right-handed 90-degree tertiary mirror), “alt-az+nasmyth-l” (“alt-az” mount with a left-handed 90-degree tertiary mirror), “phased” (antenna is “electronically steered” by summing the voltages of multiple elements, e.g. MWA), “fixed” (antenna beam pattern is fixed in azimuth and elevation, e.g., HERA), and “other” (also referred to in some formats as “bizarre”). See the “Conventions” page of the documentation for further details.

Type:

str

Parameters:
  • feed_array (array-like of str) – Feeds to define this beam for, e.g. x & y or r & l. Default is [“x”, “y”].

  • feed_angle (array-like of float) – Position angle of a given feed, units of radians. A feed angle of 0 is typically oriented toward zenith for steerable antennas, otherwise toward north for fixed antennas (e.g., HERA, LWA). More details on this can be found on the “Conventions” page of the docs.

  • mount_type (str) – Antenna mount type, which describes the optics of the antenna in question. Supported options include: “alt-az” (primary rotates in azimuth and elevation), “equatorial” (primary rotates in hour angle and declination) “orbiting” (antenna is in motion, and its orientation depends on orbital parameters), “x-y” (primary rotates first in the plane connecting east, west, and zenith, and then perpendicular to that plane), “alt-az+nasmyth-r” (“alt-az” mount with a right-handed 90-degree tertiary mirror), “alt-az+nasmyth-l” (“alt-az” mount with a left-handed 90-degree tertiary mirror), “phased” (antenna is “electronically steered” by summing the voltages of multiple elements, e.g. MWA), “fixed” (antenna beam pattern is fixed in azimuth and elevation, e.g., HERA), and “other” (also referred to in some formats as “bizarre”). See the “Conventions” page of the documentation for further details.

  • include_cross_pols (bool) – Option to include the cross polarized beams (e.g. xy and yx or en and ne) for the power beam.

  • x_orientation (str) – Physical orientation of the feed for the x feed. Not meaningful for unpolarized analytic beams, but clarifies the orientation of the dipole for for polarized beams like the ShortDipoleBeam and matches with the meaning on UVBeam objects. Used to set values of the feed_angle attribute, ignored in feed_angle is provided. Default is “east”.

yaml constructors

Analytic beams can be instantiated directly from yaml files using the !AnalyticBeam tag. The class parameter must be specified and it can be set to one of the pyuvdata provided analytic beams or to any AnalyticBeam subclass. If the subclass is not defined in pyuvdata, either the subclass must already be imported or it must be specified with the dot-pathed modules included (i.e. my_module.MyAnalyticBeam). Some analytic beams have required parameters (e.g. diameter for AiryBeams), which must also be provided, see the object definitions for details.

Some examples are provided below, note that the node key can be anything, it does not need to be beam:

A 16 meter diameter Airy beam:

beam: !AnalyticBeam
    class: AiryBeam
    diameter: 16

A classical short dipole beam (the dot-pathed module notation is not required for pyvudata beams but is shown here as an example):

beam: !AnalyticBeam
    class: pyuvdata.ShortDipoleBeam

A gaussian shaped beam with an E-Field beam sigma of 0.26 radians that has width that scales as a power law with frequency:

beam: !AnalyticBeam
    class: GaussianBeam
    reference_frequency: 120000000.
    spectral_index: -1.5
    sigma: 0.26