UVBeam
UVBeam is the main user class for primary beam models for radio telescopes. It provides import and export functionality to and from the supported file formats (beamFITS, CST, MWA primary beam) as well as methods for transforming the data (interpolating/regridding, selecting, converting types) and can be interacted with directly.
Note that there are some tricks that can help with reading in CST beam simulation files in CST Settings Files. UVBeam also provides a yaml constructor that can enable beams to be instantiated directly from yaml files (see: UVbeam yaml constructor)
Attributes
The attributes on UVBeam hold all of the metadata and data required to
describe primary beam models. Under the hood, the attributes are implemented
as properties based on pyuvdata.parameter.UVParameter
objects but
this is fairly transparent to users.
UVBeam objects can be initialized from a file using the
pyuvdata.UVBeam.from_file()
class method
(as beam = UVBeam.from_file(<filename>)
) or be initialized from arrays
in memory using the pyuvdata.UVBeam.new()
class method or as an empty
object (as beam = UVBeam()
). When an empty UVBeam object is initialized,
it has all of these attributes defined but set to None
. The attributes
can be set by reading in a data file using the pyuvdata.UVBeam.read()
method or by setting them directly on the object. Some of these attributes
are required to be set to have a fully defined data set while others are
optional. The pyuvdata.UVBeam.check()
method can be called on the
object to verify that all of the required attributes have been set in a
consistent way.
Required
These parameters are required to have a well-defined UVBeam object and are required for most kinds of beam files.
- Naxes_vec
Number of basis vectors used to represent the antenna response in each pixel. These need not align with the pixel coordinate system or even be orthogonal. The mapping of these basis vectors to directions aligned withthe pixel coordinate system is contained in the basis_vector_array. The allowed values for this parameter are 2 or 3 (or 1 if beam_type is ‘power’).
- Nfreqs
Number of frequency channels
- antenna_type
String indicating antenna type. Allowed values are “simple”, and “phased_array”
- bandpass_array
Frequency dependence of the beam. Depending on the data_normalization, this may contain only the frequency dependence of the receiving chain (‘physical’ normalization) or all the frequency dependence (‘peak’ normalization). Shape is (Nfreqs,).
- beam_type
String indicating beam type. Allowed values are ‘efield’, and ‘power’.
- data_array
Depending on beam type, either complex E-field values (‘efield’ beam type) or power values (‘power’ beam type) for beam model. Units are normalized to either peak or solid angle as given by data_normalization. The shape depends on the beam_type and pixel_coordinate_system. If it is a ‘healpix’ beam, the shape is: (Naxes_vec, Nfeeds or Npols, Nfreqs, Npixels), if it is not a healpix beam it is (Naxes_vec, Nfeeds or Npols, Nfreqs, Naxes2, Naxes1).
- data_normalization
Normalization standard of data_array, options are: “physical”, “peak” or “solid_angle”. Physical normalization means that the frequency dependence of the antenna sensitivity is included in the data_array while the frequency dependence of the receiving chain is included in the bandpass_array. Peak normalized means that for each frequency the data_arrayis separately normalized such that the peak is 1 (so the beam is dimensionless) and all direction-independent frequency dependence is moved to the bandpass_array (if the beam_type is “efield”, then peak normalized means that the absolute value of the peak is 1). Solid angle normalized means the peak normalized beam is divided by the integral of the beam over the sphere, so the beam has dimensions of 1/stradian.
- feed_name
Name of physical feed (string)
- feed_version
Version of physical feed (string)
- freq_array
Array of frequencies, center of the channel, shape (Nfreqs,), units Hz.
- history
String of history, units English
- model_name
Name of beam model (string)
- model_version
Version of beam model (string)
- pixel_coordinate_system
Pixel coordinate system, options are: “az_za”, “orthoslant_zenith”, “healpix”. “az_za” is a uniformly gridded azimuth, zenith angle coordinate system, where az runs from East to North in radians. It has axes [azimuth, zen_angle]. “orthoslant_zenith” is a orthoslant projection at zenith where y points North, x point East. It has axes [zenorth_x, zenorth_y]. “healpix” is a HEALPix map with zenith at the north pole and az, za coordinate axes (for the basis_vector_array) where az runs from East to North. It has axes [hpx_inds].
- telescope_name
Name of telescope (string)
Optional
These parameters are defined by one or more file standard but are not always required. Some of them are required depending on the beam_type, antenna_type and pixel_coordinate_systems (as noted below).
- Naxes1
Number of elements along the first pixel axis. Not required if pixel_coordinate_system is “healpix”.
- Naxes2
Number of elements along the second pixel axis. Not required if pixel_coordinate_system is “healpix”.
- Ncomponents_vec
Number of orthogonal components required to map each basis vector to vectors aligned with the pixel coordinate system. This can be equal to or smaller than Naxes_vec. The allowed values for this parameter are 2 or 3. Only required for E-field beams.
- Nelements
Required if antenna_type = “phased_array”. Number of elements in phased array
- Nfeeds
Number of feeds. Not required if beam_type is ‘power’.
- Npixels
Number of healpix pixels. Only required if pixel_coordinate_system is ‘healpix’.
- Npols
Number of polarizations. Only required if beam_type is “power”.
- axis1_array
Coordinates along first pixel axis. Not required if pixel_coordinate_system is “healpix”.
- axis2_array
Coordinates along second pixel axis. Not required if pixel_coordinate_system is “healpix”.
- basis_vector_array
Beam basis vector components, essentially the mapping between the directions that the electrical field values are recorded in to the directions aligned with the pixel coordinate system (or azimuth/zenith angle for HEALPix beams).Not required if beam_type is “power”. The shape depends on the pixel_coordinate_system, if it is “healpix”, the shape is: (Naxes_vec, Ncomponents_vec, Npixels), otherwise it is (Naxes_vec, Ncomponents_vec, Naxes2, Naxes1)
- coupling_matrix
Required if antenna_type = “phased_array”. Matrix of complex element couplings, units: dB, shape: (Nelements, Nelements, Nfeeds, Nfeeds, Nfreqs).
- delay_array
Required if antenna_type = “phased_array”. Array of element delays, units: seconds, shape: (Nelements)
- element_coordinate_system
Required if antenna_type = “phased_array”. Element coordinate system, options are: n-e or x-y
- element_location_array
Required if antenna_type = “phased_array”. Array of element locations in element coordinate system, shape: (2, Nelements)
- extra_keywords
Any user supplied extra keywords, type=dict. Keys should be 8 character or less strings if writing to beam fits files. Use the special key “comment” for long multi-line string comments.
- feed_array
Array of feed orientations. shape (Nfeeds). options are: n/e or x/y or r/l. Not required if beam_type is “power”.
- filename
List of strings containing the unique basenames (not the full path) of input files.
- gain_array
Required if antenna_type = “phased_array”. Array of element gains, units: dB, shape: (Nelements)
- loss_array
Array of antenna losses, units dB? Shape (Nfreqs,).
- mismatch_array
Array of antenna-amplifier mismatches, units ? Shape (Nfreqs,).
- nside
Healpix nside parameter. Only required if pixel_coordinate_system is ‘healpix’.
- ordering
Healpix ordering parameter, allowed values are “ring” and “nested”. Only required if pixel_coordinate_system is “healpix”.
- pixel_array
Healpix pixel numbers. Only required if pixel_coordinate_system is ‘healpix’.
- polarization_array
Array of polarization integers, shape (Npols). Uses the same convention as UVData: pseudo-stokes 1:4 (pI, pQ, pU, pV); circular -1:-4 (RR, LL, RL, LR); linear -5:-8 (XX, YY, XY, YX). Only required if beam_type is “power”.
- receiver_temperature_array
Array of receiver temperatures, units K. Shape (Nfreqs,).
- reference_impedance
Reference impedance of the beam model. The radiated E-farfield or the realised gain depend on the impedance of the port used to excite the simulation. This is the reference impedance (Z0) of the simulation. units: Ohms
- s_parameters
S parameters of receiving chain, ordering: s11, s12, s21, s22. see https://en.wikipedia.org/wiki/Scattering_parameters#Two-Port_S-ParametersShape (4, Nfreqs).
- x_orientation
Orientation of the physical dipole corresponding to what is labelled as the x polarization. Options are “east” (indicating east/west orientation) and “north” (indicating north/south orientation)
Methods
- class pyuvdata.UVBeam[source]
A class for defining a radio telescope antenna beam.
- UVParameter objects
For full list see the documentation on ReadTheDocs: http://pyuvdata.readthedocs.io/en/latest/. Some are always required, some are required for certain beam_types, antenna_types and pixel_coordinate_systems and others are always optional.
- static new(**kwargs)[source]
Create a new UVBeam object.
All parameters are passed through to the
new_uvbeam()
function.- Returns:
UVBeam – A new UVBeam object.
- Other Parameters:
telescope_name (str) – Telescope name.
data_normalization (str) – Normalization standard of data_array, options are: “physical”, “peak” or “solid_angle”. Physical normalization means that the frequency dependence of the antenna sensitivity is included in the data_array while the frequency dependence of the receiving chain is included in the bandpass_array. Peak normalized means that for each frequency the data_array is separately normalized such that the peak is 1 (so the beam is dimensionless) and all direction-independent frequency dependence is moved to the bandpass_array (if the beam_type is “efield”, then peak normalized means that the absolute value of the peak is 1). Solid angle normalized means the peak normalize beam is divided by the integral of the beam over the sphere, so the beam has dimensions of 1/steradian.
freq_array (ndarray of float) – Array of frequencies in Hz.
feed_name (str) – Name of the physical feed.
feed_version (str) – Version of the physical feed.
model_name (str) – Name of the beam model.
model_version (str) – Version of the beam model.
feed_array (ndarray of str) – Array of feed orientations. Options are: n/e or x/y or r/l. Must be provided for an E-field beam.
polarization_array (ndarray of str or int) – Array of polarization integers or strings (eg. ‘xx’ or ‘ee’). Must be provided for a power beam.
x_orientation (str, optional) – Orientation of the x-axis. Options are ‘east’, ‘north’, ‘e’, ‘n’, ‘ew’, ‘ns’.
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.axis1_array (ndarray of float) – Coordinates along first pixel axis (e.g. azimuth for an azimuth/zenith angle coordinate system). 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). 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.
basis_vector_array (ndarray of float) – Beam basis vector components, essentially the mapping between the directions that the electrical field values are recorded in to the directions aligned with the pixel coordinate system (or azimuth/zenith angle for HEALPix beams). Defaults to unit vectors aligned with the pixel coordinate systems for E-field beams.
bandpass_array (ndarray of float) – Frequency dependence of the beam. Depending on the
data_normalization
this may contain only the frequency dependence of the receiving chain (“physical” normalization) or all the frequency dependence (“peak” normalization). Must be the same length as thefreq_array
. Defaults to an array of all ones the length of thefreq_array
.element_location_array (ndarray of float) – Array of phased array element locations in the element_coordinate_system. Must be a 2 dimensional array where the first dimension indexes the coordinate system (so should be length 2) and the second dimension indexes the phased array elements. Only used for phase array antennas.
element_coordinate_system (str) – Coordinate system for describing the layout of a phased array feed. Options are: n-e or x-y. Defaults to “x-y” if
element_location_array
is provided. Only used for phase array antennas.delay_array (ndarray of float) – Array of element delays in seconds. Defaults to an array of all zeros if
element_location_array
is provided. Only used for phase array antennas.gain_array (ndarray of float) – Array of element gains in dB. Defaults to an array of all ones if
element_location_array
is provided. Only used for phase array antennas.coupling_matrix (ndarray of float) – Matrix of complex element couplings in dB. Must be an ndarray of shape (Nelements, Nelements, Nfeeds, Nfeeds, Nfreqs). Defaults to an array with couplings of one for the same element and zero for other elements if
element_location_array
is provided. Only used for phase array antennas.data_array (ndarray of float, optional) – Either complex E-field values (if feed_array is given) or power values (if polarization_array is given) for the beam model. Units are normalized to either peak or solid angle as given by data_normalization. If None (the default), the data_array is initialized to an array of the appropriate shape and type containing all zeros.
history (str, optional) – History string to be added to the object. Default is a simple string containing the date and time and pyuvdata version.
**kwargs – All other keyword arguments are added to the object as attributes.
- check(*, check_extra=True, run_check_acceptability=True, check_auto_power=False, fix_auto_power=False)[source]
Check that all required parameters are set reasonably.
Check that required parameters exist and have appropriate shapes. Optionally check if the values are acceptable.
- Parameters:
check_extra (bool) – Option to check optional parameters as well as required ones.
run_check_acceptability (bool) – Option to check if values in required parameters are acceptable.
check_auto_power (bool) – For power beams, check whether the auto polarization beams have non-zero imaginary values in the data_array (which should not mathematically exist).
fix_auto_power (bool) – For power beams, if auto polarization beams with imaginary values are found, fix those values so that they are real-only in data_array. Ignored if check_auto_power is False.
- efield_to_power(*, calc_cross_pols=True, keep_basis_vector=False, run_check=True, check_extra=True, run_check_acceptability=True, inplace=True)[source]
Convert E-field beam to power beam.
- Parameters:
calc_cross_pols (bool) – If True, calculate the crossed polarization beams (e.g. ‘xy’ and ‘yx’), otherwise only calculate the same polarization beams (e.g. ‘xx’ and ‘yy’).
keep_basis_vector (bool) – If True, keep the directionality information and just multiply the efields for each basis vector separately (caution: this is not what is standardly meant by the power beam).
inplace (bool) – Option to apply conversion directly on self or to return a new UVBeam object.
run_check (bool) – Option to check for the existence and proper shapes of the required parameters after converting to power.
run_check_acceptability (bool) – Option to check acceptable range of the values of required parameters after converting to power.
check_extra (bool) – Option to check optional parameters as well as required ones.
- efield_to_pstokes(*, inplace=True, run_check=True, check_extra=True, run_check_acceptability=True)[source]
Convert E-field to pseudo-stokes power.
Following https://arxiv.org/pdf/1802.04151.pdf, using the equation:
M_ij = Tr(sigma_i J sigma_j J^*)
where sigma_i and sigma_j are Pauli matrices.
- Parameters:
inplace (bool) – Option to apply conversion directly on self or to return a new UVBeam object.
run_check (bool) – Option to check for the existence and proper shapes of the required parameters after converting to power.
run_check_acceptability (bool) – Option to check acceptable range of the values of required parameters after converting to power.
check_extra (bool) – Option to check optional parameters as well as required ones.
- interp(*, az_array=None, za_array=None, interpolation_function=None, freq_interp_kind=None, az_za_grid=False, healpix_nside=None, healpix_inds=None, freq_array=None, freq_interp_tol=1.0, polarizations=None, return_bandpass=False, return_coupling=False, return_basis_vector: bool | None = None, reuse_spline=False, spline_opts=None, new_object=False, run_check=True, check_extra=True, run_check_acceptability=True, check_azza_domain: bool = True)[source]
Interpolate beam to given frequency, az & za locations or Healpix pixel centers.
Uses the function specified in interpolation_function, which defaults to “az_za_simple” for objects with the “az_za” pixel_coordinate_system and “healpix_simple” for objects with the “healpix” pixel_coordinate_system. Currently supported interpolation functions include:
“az_za_simple”: Uses scipy RectBivariate spline interpolation, can only be used on objects with an “az_za” pixel_coordinate_system.
“az_za_map_coordinates”: Uses scipy map_coordinates interpolation, can only be used on objects with an “az_za” pixel_coordinate_system.
“healpix_simple”: Uses HEALPix nearest-neighbor bilinear interpolation, can only be used on objects with a “healpix” pixel_coordinate_system.
- Parameters:
az_array (array_like of floats, optional) – Azimuth values to interpolate to in radians, either specifying the azimuth positions for every interpolation point or specifying the azimuth vector for a meshgrid if az_za_grid is True.
za_array (array_like of floats, optional) – Zenith values to interpolate to in radians, either specifying the zenith positions for every interpolation point or specifying the zenith vector for a meshgrid if az_za_grid is True.
interpolation_function (str, optional) – Specify the interpolation function to use. Defaults to: “az_za_simple” for objects with the “az_za” pixel_coordinate_system and “healpix_simple” for objects with the “healpix” pixel_coordinate_system. “az_za_map_coordinates” is also available for objects with the “az_za” pixel_coordinate_system.
freq_interp_kind (str) – Interpolation method to use frequency. See scipy.interpolate.interp1d for details. Defaults to “cubic” (Note that this is a change. It used to default to “linear” when it was assigned to the object. However, multiple groups have found that a linear interpolation leads to nasty artifacts in visibility simulations for EoR applications.)
az_za_grid (bool) – Option to treat the az_array and za_array as the input vectors for points on a mesh grid.
healpix_nside (int, optional) – HEALPix nside parameter if interpolating to HEALPix pixels.
healpix_inds (array_like of int, optional) – HEALPix indices to interpolate to. Defaults to all indices in the map if healpix_nside is set and az_array and za_array are None.
freq_array (array_like of floats, optional) – Frequency values to interpolate to.
freq_interp_tol (float) – Frequency distance tolerance [Hz] of nearest neighbors. If all elements in freq_array have nearest neighbor distances within the specified tolerance then return the beam at each nearest neighbor, otherwise interpolate the beam.
polarizations (list of str) – polarizations to interpolate if beam_type is ‘power’. Default is all polarizations in self.polarization_array.
return_bandpass (bool) – Option to return the bandpass. Only applies if new_object is False.
return_coupling (bool) – Option to return the interpolated coupling matrix, only applies if antenna_type is “phased_array” and new_object is False.
new_object (bool) – Option to return a new UVBeam object with the interpolated data, if possible. Note that this is only possible for Healpix pixels or if az_za_grid is True and az_array and za_array are evenly spaced or for frequency only interpolation.
reuse_spline (bool) – Save the interpolation functions for reuse. Only applies for az_za_simple and az_za_map_coordinates interpolation.
spline_opts (dict) – Provide options to numpy.RectBivariateSpline. This includes spline order parameters kx and ky, and smoothing parameter s. Applies for az_za_simple and az_za_map_coordinates interpolation.
run_check (bool) – Only used if new_object is True. Option to check for the existence and proper shapes of required parameters on the new object.
check_extra (bool) – Only used if new_object is True. Option to check optional parameters as well as required ones on the new object.
run_check_acceptability (bool) – Only used if new_object is True. Option to check acceptable range of the values of required parameters on the new object.
check_azza_domain (bool) – Whether to check the domain of az/za to ensure that they are covered by the intrinsic data array. Checking them can be quite computationally expensive. Conversely, if the passed az/za are outside of the domain, they will be silently extrapolated and the behavior is not well-defined. Only applies for az_za_simple and az_za_map_coordinates interpolation.
return_basis_vector (bool) – Whether to return the interpolated basis vectors. Prior to v3.1.1 these were always returned. In v3.3+ they will _not_ be returned by default (and not computed by default, unless new_object=True).
- Returns:
array_like of float or complex or a UVBeam object – Either an array of interpolated values or a UVBeam object if new_object is True. The shape of the interpolated data will be: (Naxes_vec, Nfeeds or Npols, Nfreqs or freq_array.size if freq_array is passed, Npixels/(Naxis1, Naxis2) or az_array.size if az/za_arrays are passed) or (Naxes_vec, Nfeeds or Npols, Nfreqs or freq_array.size if freq_array is passed, Npixels/(Naxis1, Naxis2)
interp_basis_vector (array_like of float, optional) – The array of interpolated basis vectors (or self.basis_vector_array if az/za_arrays are not passed). Only returned if new_object is False. shape: (Naxes_vec, Ncomponents_vec, Npixels/(Naxis1, Naxis2) or az_array.size if az/za_arrays are passed)
interp_bandpass (array_like of float, optional) – The interpolated bandpass, only returned if return_bandpass is True and new_object is False. Shape: (freq_array.size,)
interp_coupling_matrix (array_like of complex, optional) – The interpolated coupling matrix. Only returned if return_coupling is True and new_object is False. Shape: (Nelements, Nelements, Nfeeds, Nfeeds, freq_array.size)
- to_healpix(*, nside=None, interpolation_function=None, run_check=True, check_extra=True, run_check_acceptability=True, inplace=True)[source]
Convert beam to the healpix coordinate system using interpolation.
Note that this interpolation isn’t perfect. Interpolating an Efield beam and then converting to power gives a different result than converting to power and then interpolating at about a 5% level.
Uses the function specified in interpolation_function, defaults to the “az_za_simple” for objects with the “az_za” pixel_coordinate_system. Currently supported interpolation functions for beams that are not already in a healpix coordinate system include:
“az_za_simple”: Uses scipy RectBivariate spline interpolation, can only be used on objects with an “az_za” pixel_coordinate_system.
“az_za_map_coordinates”: Uses scipy map_coordinates interpolation, can only be used on objects with an “az_za” pixel_coordinate
- Parameters:
nside (int) – The nside to use for the Healpix map. If not specified, use the nside that gives the closest resolution that is higher than the input resolution.
interpolation_function (str, optional) – Specify the interpolation function to use. Defaults to to: “az_za_simple” for objects with the “az_za” pixel_coordinate_system. “az_za_map_coordinates” is also available for objects with the “az_za” pixel_coordinate_system.
run_check (bool) – Option to check for the existence and proper shapes of required parameters after converting to healpix.
check_extra (bool) – Option to check optional parameters as well as required ones.
run_check_acceptability (bool) – Option to check acceptable range of the values of required parameters after combining objects
inplace (bool) – Option to perform the interpolation directly on self or return a new UVBeam object.
- get_beam_area(pol='pI')[source]
Compute the integral of the beam in units of steradians.
Pseudo-Stokes ‘pI’ (I), ‘pQ’(Q), ‘pU’(U), ‘pV’(V) beam and linear dipole ‘XX’, ‘XY’, ‘YX’ and ‘YY’ are supported. See Equations 4 and 5 of Moore et al. (2017) ApJ 836, 154 or arxiv:1502.05072 and Kohn et al. (2018) or https://arxiv.org/pdf/1802.04151.pdf for details.
- Parameters:
pol (str or int) – polarization string or integer, Ex. a pseudo-stokes pol ‘pI’, or a linear pol ‘XX’.
- Returns:
omega (float) – Integral of the beam across the sky, units: steradians.
- get_beam_sq_area(pol='pI')[source]
Compute the integral of the beam^2 in units of steradians.
Pseudo-Stokes ‘pI’ (I), ‘pQ’(Q), ‘pU’(U), ‘pV’(V) beam and linear dipole ‘XX’, ‘XY’, ‘YX’ and ‘YY’ are supported. See Equations 4 and 5 of Moore et al. (2017) ApJ 836, 154 or arxiv:1502.05072 for details.
- Parameters:
pol (str or int) – polarization string or integer, Ex. a pseudo-stokes pol ‘pI’, or a linear pol ‘XX’.
- Returns:
omega (float) – Integral of the beam^2 across the sky, units: steradians.
- select(*, axis1_inds=None, axis2_inds=None, pixels=None, frequencies=None, freq_chans=None, feeds=None, polarizations=None, inplace=True, run_check=True, check_extra=True, run_check_acceptability=True)[source]
Downselect data to keep on the object along various axes.
Axes that can be selected along include image axis indices or pixels (if healpix), frequencies and feeds or polarizations (if power).
The history attribute on the object will be updated to identify the operations performed.
- Parameters:
axis1_indss (array_like of int, optional) – The indices along the first image axis to keep in the object. Cannot be set if pixel_coordinate_system is “healpix”.
axis2_inds (array_like of int, optional) – The indices along the second image axis to keep in the object. Cannot be set if pixel_coordinate_system is “healpix”.
pixels (array_like of int, optional) – The healpix pixels to keep in the object. Cannot be set if pixel_coordinate_system is not “healpix”.
frequencies (array_like of float, optional) – The frequencies to keep in the object.
freq_chans (array_like of int, optional) – The frequency channel numbers to keep in the object.
feeds (array_like of str, optional) – The feeds to keep in the object. If the x_orientation attribute is set, the physical dipole strings (e.g. “n”, “e”) are also supported. Cannot be set if the beam_type is “power”.
polarizations (array_like of int or str, optional) – The polarizations to keep in the object. Cannot be set if the beam_type is “efield”. If passing strings, the canonical polarization strings (e.g. “xx”, “rr”) are supported and if the x_orientation attribute is set, the physical dipole strings (e.g. “nn”, “ee”) are also supported.
inplace (bool) – Option to perform the select directly on self or return a new UVBeam object, which is a subselection of self.
run_check (bool) – Option to check for the existence and proper shapes of required parameters after downselecting data on this object.
check_extra (bool) – Option to check optional parameters as well as required ones.
run_check_acceptability (bool) – Option to check acceptable range of the values of required parameters after downselecting data on this object.
- read_beamfits(filename, **kwargs)[source]
Read in data from a beamfits file.
- Parameters:
filename (str or list of str) – The beamfits file or list of files to read from.
use_future_array_shapes (bool) – Defunct option, will result in an error in version 3.2.
run_check (bool) – Option to check for the existence and proper shapes of required parameters after reading in the file.
check_extra (bool) – Option to check optional parameters as well as required ones.
run_check_acceptabilit (bool) – Option to check acceptable range of the values of required parameters after reading in the file.
check_auto_power (bool) – For power beams, check whether the auto polarization beams have non-zero imaginary values in the data_array (which should not mathematically exist).
fix_auto_power (bool) – For power beams, if auto polarization beams with imaginary values are found, fix those values so that they are real-only in data_array.
freq_range (tuple of float in Hz) – If given, the lower and upper limit of the frequencies to read in. Default is to read in all frequencies. Restricting the frequencies reduces peak memory usage.
az_range (tuple of float in deg) – The azimuth range to read in, if the beam is specified in az/za coordinates. Default is to read in all azimuths. Restricting the azimuth reduces peak memory usage.
za_range (tuple of float in deg) – The zenith angle range to read in, if the beam is specified in za/za coordinates. Default is to read in all za. Restricting the za reduces peak memory.
- read_cst_beam(filename, *, beam_type='power', use_future_array_shapes=None, feed_pol=None, rotate_pol=None, frequency=None, telescope_name=None, feed_name=None, feed_version=None, model_name=None, model_version=None, history=None, x_orientation=None, reference_impedance=None, extra_keywords=None, frequency_select=None, run_check=True, check_extra=True, run_check_acceptability=True, check_auto_power=True, fix_auto_power=True)[source]
Read in data from a cst file.
- Parameters:
filename (str) – Either a settings yaml file or a cst text file or list of cst text files to read from. If a list is passed, the files are combined along the appropriate axes.
Settings yaml files must include the following keywords:
- telescope_name (str)- feed_name (str)- feed_version (str)- model_name (str)- model_version (str)- history (str)- frequencies (list(float))- cst text filenames (list(str)) – path relative to yaml file- feed_pol (str) or (list(str))and they may include the following optional keywords:
- x_orientation (str): Optional but strongly encouraged!- ref_imp (float): beam model reference impedance- sim_beam_type (str): e.g. ‘E-farfield’- all other fields will go into the extra_keywords attributeMore details and an example are available in the docs (cst_settings_yaml.rst). Specifying any of the associated keywords to this function will override the values in the settings file.
beam_type (str) – What beam_type to read in (‘power’ or ‘efield’).
use_future_array_shapes (bool) – Defunct option, will result in an error in version 3.2.
feed_pol (str) – The feed or polarization or list of feeds or polarizations the files correspond to. Defaults to ‘x’ (meaning x for efield or xx for power beams).
rotate_pol (bool) – If True, assume the structure in the simulation is symmetric under 90 degree rotations about the z-axis (so that the y polarization can be constructed by rotating the x polarization or vice versa). Default: True if feed_pol is a single value or a list with all the same values in it, False if it is a list with varying values.
frequency (float or list of float, optional) – The frequency or list of frequencies corresponding to the filename(s). This is assumed to be in the same order as the files. If not passed, the code attempts to parse it from the filenames.
telescope_name (str, optional) – The name of the telescope corresponding to the filename(s).
feed_name (str, optional) – The name of the feed corresponding to the filename(s).
feed_version (str, optional) – The version of the feed corresponding to the filename(s).
model_name (str, optional) – The name of the model corresponding to the filename(s).
model_version (str, optional) – The version of the model corresponding to the filename(s).
history (str, optional) – A string detailing the history of the filename(s).
x_orientation (str, optional) – Orientation of the physical dipole corresponding to what is labelled as the x polarization. Options are “east” (indicating east/west orientation) and “north” (indicating north/south orientation)
reference_impedance (float, optional) – The reference impedance of the model(s).
extra_keywords (dict, optional) – A dictionary containing any extra_keywords.
frequency_select (list of float, optional) – Only used if the file is a yaml file. Indicates which frequencies to include (only read in files for those frequencies)
run_check (bool) – Option to check for the existence and proper shapes of required parameters after reading in the file.
check_extra (bool) – Option to check optional parameters as well as required ones.
run_check_acceptability (bool) – Option to check acceptable range of the values of required parameters after reading in the file.
check_auto_power (bool) – For power beams, check whether the auto polarization beams have non-zero imaginary values in the data_array (which should not mathematically exist).
fix_auto_power (bool) – For power beams, if auto polarization beams with imaginary values are found, fix those values so that they are real-only in data_array.
- read_mwa_beam(h5filepath, **kwargs)[source]
Read in the full embedded element MWA beam.
Note that the azimuth convention in for the UVBeam object is different than the azimuth convention in the mwa_pb repo. In that repo, the azimuth convention is changed from the native FEKO convention (the FEKO convention is the same as the UVBeam convention). The convention in the mwa_pb repo has a different zero point and a different direction (so it is in a left handed coordinate system).
- Parameters:
h5filepath (str) – path to input h5 file containing the MWA full embedded element spherical harmonic modes. Download via wget http://ws.mwatelescope.org/static/mwa_full_embedded_element_pattern.h5 (This reader is based on https://github.com/MWATelescope/mwa_pb).
use_future_array_shapes (bool) – Defunct option, will result in an error in version 3.2.
delays (array of ints) – Array of MWA beamformer delay steps. Should be shape (n_pols, n_dipoles).
amplitudes (array of floats) – Array of dipole amplitudes, these are absolute values (i.e. relatable to physical units). Should be shape (n_pols, n_dipoles).
pixels_per_deg (float) – Number of theta/phi pixels per degree. Sets the resolution of the beam.
freq_range (array_like of float) – Range of frequencies to include in Hz, defaults to all available frequencies. Must be length 2.
run_check (bool) – Option to check for the existence and proper shapes of required parameters after reading in the file.
check_extra (bool) – Option to check optional parameters as well as required ones.
run_check_acceptability (bool) – Option to check acceptable range of the values of required parameters after reading in the file.
check_auto_power (bool) – For power beams, check whether the auto polarization beams have non-zero imaginary values in the data_array (which should not mathematically exist).
fix_auto_power (bool) – For power beams, if auto polarization beams with imaginary values are found, fix those values so that they are real-only in data_array.
- read(filename, *, file_type=None, skip_bad_files=False, use_future_array_shapes=None, run_check=True, check_extra=True, run_check_acceptability=True, check_auto_power=True, fix_auto_power=True, az_range=None, za_range=None, freq_range=None, beam_type='power', feed_pol=None, rotate_pol=None, frequency=None, telescope_name=None, feed_name=None, feed_version=None, model_name=None, model_version=None, history=None, x_orientation=None, reference_impedance=None, extra_keywords=None, frequency_select=None, delays=None, amplitudes=None, pixels_per_deg=5)[source]
Read a generic file into a UVBeam object.
Some parameters only apply to certain file types.
- Parameters:
filename (str or array_like of str) – The file(s) or list(s) (or array(s)) of files to read from.
For cst yaml files only:
Settings yaml files must include the following keywords:
- telescope_name (str)- feed_name (str)- feed_version (str)- model_name (str)- model_version (str)- history (str)- frequencies (list(float))- cst text filenames (list(str)) – path relative to yaml file- feed_pol (str) or (list(str))and they may include the following optional keywords:
- x_orientation (str): Optional but strongly encouraged!- ref_imp (float): beam model reference impedance- sim_beam_type (str): e.g. ‘E-farfield’- all other fields will go into the extra_keywords attributeMore details and an example are available in the docs (cst_settings_yaml.rst). Specifying any of the associated keywords to this function will override the values in the settings file.
file_type (str) – One of [‘mwa_beam’, ‘beamfits’, ‘cst’] or None. If None, the code attempts to guess what the file type is. based on file extensions (mwa_beam: .hdf5, .h5; cst: .yaml, .txt; beamfits: .fits, .beamfits). Note that if a list of datasets is passed, the file type is determined from the first dataset.
skip_bad_files (bool) – Option when reading multiple files to catch read errors such that the read continues even if one or more files are corrupted. Files that produce errors will be printed. Default is False (files will not be skipped).
use_future_array_shapes (bool) – Defunct option, will result in an error in version 3.2.
- Checking:
run_check (bool) – Option to check for the existence and proper shapes of parameters after after reading in the file (the default is True, meaning the check will be run). Ignored if read_data is False.
check_extra (bool) – Option to check optional parameters as well as required ones (the default is True, meaning the optional parameters will be checked). Ignored if read_data is False.
run_check_acceptability (bool) – Option to check acceptable range of the values of parameters after reading in the file (the default is True, meaning the acceptable range check will be done). Ignored if read_data is False.
check_auto_power (bool) – For power beams, check whether the auto polarization beams have non-zero imaginary values in the data_array (which should not mathematically exist).
fix_auto_power (bool) – For power beams, if auto polarization beams with imaginary values are found, fix those values so that they are real-only in data_array.
Beamfits
——–
az_range (tuple of float in deg) – The azimuth range to read in, if the beam is specified in az/za coordinates. Default is to read in all azimuths. Restricting the azimuth reduces peak memory usage. Only used for beamfits files that have their coordinates in az/za grid.
za_range (tuple of float in deg) – The zenith angle range to read in, if the beam is specified in za/za coordinates. Default is to read in all za. Restricting the za reduces peak memory. Only used for beamfits files that have their coordinates in az/za grid.
Beamfits & MWA
————–
freq_range (array_like of float) – Range of frequencies to include in Hz, defaults to all available frequencies. Must be length 2. Only applies to mwa_beam and beamfits type files. For beamfits, this will cause a partial read (i.e. reduce peak memory usage).
CST
—
beam_type (str) – What beam_type to read in (‘power’ or ‘efield’). Only applies to cst file types.
feed_pol (str) – The feed or polarization or list of feeds or polarizations the files correspond to. Defaults to ‘x’ (meaning x for efield or xx for power beams). Only applies to cst file types.
rotate_pol (bool) – If True, assume the structure in the simulation is symmetric under 90 degree rotations about the z-axis (so that the y polarization can be constructed by rotating the x polarization or vice versa). Default: True if feed_pol is a single value or a list with all the same values in it, False if it is a list with varying values. Only applies to cst file types.
frequency (float or list of float, optional) – The frequency or list of frequencies corresponding to the filename(s). This is assumed to be in the same order as the files. If not passed, the code attempts to parse it from the filenames. Only applies to cst file types.
telescope_name (str, optional) – The name of the telescope corresponding to the filename(s). Only applies to cst file types.
feed_name (str, optional) – The name of the feed corresponding to the filename(s). Only applies to cst file types.
feed_version (str, optional) – The version of the feed corresponding to the filename(s). Only applies to cst file types.
model_name (str, optional) – The name of the model corresponding to the filename(s). Only applies to cst file types.
model_version (str, optional) – The version of the model corresponding to the filename(s). Only applies to cst file types.
history (str, optional) – A string detailing the history of the filename(s). Only applies to cst file types.
x_orientation (str, optional) – Orientation of the physical dipole corresponding to what is labelled as the x polarization. Options are “east” (indicating east/west orientation) and “north” (indicating north/south orientation) Only applies to cst file types.
reference_impedance (float, optional) – The reference impedance of the model(s). Only applies to cst file types.
extra_keywords (dict, optional) – A dictionary containing any extra_keywords. Only applies to cst file types.
frequency_select (list of float, optional) – Only used if the file is a yaml file. Indicates which frequencies to include (only read in files for those frequencies) Only applies to cst file types.
MWA
—
delays (array of ints) – Array of MWA beamformer delay steps. Should be shape (n_pols, n_dipoles). Only applies to mwa_beam type files.
amplitudes (array of floats) – Array of dipole amplitudes, these are absolute values (i.e. relatable to physical units). Should be shape (n_pols, n_dipoles). Only applies to mwa_beam type files.
pixels_per_deg (float) – Number of theta/phi pixels per degree. Sets the resolution of the beam. Only applies to mwa_beam type files.
- Raises:
ValueError – If the file_type is not set and cannot be determined from the file name.
- classmethod from_file(filename, **kwargs)[source]
Initialize a new UVBeam object by reading the input file(s).
Some parameters only apply to certain file types.
- Parameters:
filename (str or array_like of str) – The file(s) or list(s) (or array(s)) of files to read from.
For cst yaml files only:
Settings yaml files must include the following keywords:
- telescope_name (str)- feed_name (str)- feed_version (str)- model_name (str)- model_version (str)- history (str)- frequencies (list(float))- cst text filenames (list(str)) – path relative to yaml file- feed_pol (str) or (list(str))and they may include the following optional keywords:
- x_orientation (str): Optional but strongly encouraged!- ref_imp (float): beam model reference impedance- sim_beam_type (str): e.g. ‘E-farfield’- all other fields will go into the extra_keywords attributeMore details and an example are available in the docs (cst_settings_yaml.rst). Specifying any of the associated keywords to this function will override the values in the settings file.
file_type (str) – One of [‘mwa_beam’, ‘beamfits’, ‘cst’] or None. If None, the code attempts to guess what the file type is. based on file extensions (mwa_beam: .hdf5, .h5; cst: .yaml, .txt; beamfits: .fits, .beamfits). Note that if a list of datasets is passed, the file type is determined from the first dataset.
skip_bad_files (bool) – Option when reading multiple files to catch read errors such that the read continues even if one or more files are corrupted. Files that produce errors will be printed. Default is False (files will not be skipped).
use_future_array_shapes (bool) – Defunct option, will result in an error in version 3.2.
- Checking:
run_check (bool) – Option to check for the existence and proper shapes of parameters after after reading in the file (the default is True, meaning the check will be run). Ignored if read_data is False.
check_extra (bool) – Option to check optional parameters as well as required ones (the default is True, meaning the optional parameters will be checked). Ignored if read_data is False.
run_check_acceptability (bool) – Option to check acceptable range of the values of parameters after reading in the file (the default is True, meaning the acceptable range check will be done). Ignored if read_data is False.
check_auto_power (bool) – For power beams, check whether the auto polarization beams have non-zero imaginary values in the data_array (which should not mathematically exist).
fix_auto_power (bool) – For power beams, if auto polarization beams with imaginary values are found, fix those values so that they are real-only in data_array.
Beamfits
——–
az_range (tuple of float in deg) – The azimuth range to read in, if the beam is specified in az/za coordinates. Default is to read in all azimuths. Restricting the azimuth reduces peak memory usage. Only used for beamfits files that have their coordinates in az/za grid.
za_range (tuple of float in deg) – The zenith angle range to read in, if the beam is specified in za/za coordinates. Default is to read in all za. Restricting the za reduces peak memory. Only used for beamfits files that have their coordinates in az/za grid.
Beamfits & MWA
————–
freq_range (array_like of float) – Range of frequencies to include in Hz, defaults to all available frequencies. Must be length 2. Only applies to mwa_beam and beamfits type files. For beamfits, this will cause a partial read (i.e. reduce peak memory usage).
CST
—
beam_type (str) – What beam_type to read in (‘power’ or ‘efield’). Only applies to cst file types.
feed_pol (str) – The feed or polarization or list of feeds or polarizations the files correspond to. Defaults to ‘x’ (meaning x for efield or xx for power beams). Only applies to cst file types.
rotate_pol (bool) – If True, assume the structure in the simulation is symmetric under 90 degree rotations about the z-axis (so that the y polarization can be constructed by rotating the x polarization or vice versa). Default: True if feed_pol is a single value or a list with all the same values in it, False if it is a list with varying values. Only applies to cst file types.
frequency (float or list of float, optional) – The frequency or list of frequencies corresponding to the filename(s). This is assumed to be in the same order as the files. If not passed, the code attempts to parse it from the filenames. Only applies to cst file types.
telescope_name (str, optional) – The name of the telescope corresponding to the filename(s). Only applies to cst file types.
feed_name (str, optional) – The name of the feed corresponding to the filename(s). Only applies to cst file types.
feed_version (str, optional) – The version of the feed corresponding to the filename(s). Only applies to cst file types.
model_name (str, optional) – The name of the model corresponding to the filename(s). Only applies to cst file types.
model_version (str, optional) – The version of the model corresponding to the filename(s). Only applies to cst file types.
history (str, optional) – A string detailing the history of the filename(s). Only applies to cst file types.
x_orientation (str, optional) – Orientation of the physical dipole corresponding to what is labelled as the x polarization. Options are “east” (indicating east/west orientation) and “north” (indicating north/south orientation) Only applies to cst file types.
reference_impedance (float, optional) – The reference impedance of the model(s). Only applies to cst file types.
extra_keywords (dict, optional) – A dictionary containing any extra_keywords. Only applies to cst file types.
frequency_select (list of float, optional) – Only used if the file is a yaml file. Indicates which frequencies to include (only read in files for those frequencies) Only applies to cst file types.
MWA
—
delays (array of ints) – Array of MWA beamformer delay steps. Should be shape (n_pols, n_dipoles). Only applies to mwa_beam type files.
amplitudes (array of floats) – Array of dipole amplitudes, these are absolute values (i.e. relatable to physical units). Should be shape (n_pols, n_dipoles). Only applies to mwa_beam type files.
pixels_per_deg (float) – Number of theta/phi pixels per degree. Sets the resolution of the beam. Only applies to mwa_beam type files.
- Raises:
ValueError – If the file_type is not set and cannot be determined from the file name.
- write_beamfits(filename, **kwargs)[source]
Write the data to a beamfits file.
- Parameters:
filename (str) – The beamfits file to write to.
run_check (bool) – Option to check for the existence and proper shapes of required parameters before writing the file.
check_extra (bool) – Option to check optional parameters as well as required ones.
run_check_acceptability (bool) – Option to check acceptable range of the values of required parameters before writing the file.
check_auto_power (bool) – For power beams, check whether the auto polarization beams have non-zero imaginary values in the data_array (which should not mathematically exist).
fix_auto_power (bool) – For power beams, if auto polarization beams with imaginary values are found, fix those values so that they are real-only in data_array.
clobber (bool) – Option to overwrite the filename if the file already exists.
UVBeam yaml constructor
UVbeams can be instantiated directly from yaml files using the
!UVBeam
tag. The filename
parameter must be specified and
any other parameter that can be passed to the pyuvdata.UVBeam.read()
method can also be specified.
Some examples are provided below, note that the node key can be anything,
it does not need to be beam
:
A simple UVBeam specification:
beam: !UVBeam
filename: hera.beamfits
An UVbeam specification with some keywords to pass to UVBeam.read
:
beam: !UVBeam
filename: mwa_full_EE_test.h5
pixels_per_deg: 1
freq_range: [100.e+6, 200.e+6]
CST Settings Files
The text files saved out of CST beam simulations do not have much of the critical metadata needed for UVBeam objects. This required metadata can be set via keywords when the files are read in, but it is better for the metadata to be specified once and carried with the data files. To that end, we developed a yaml settings file specification to carry all the metadata. This format is very human readable and writeable and we encourage using such a file as the best way to ensure the metadata is preserved.
Required Fields
The following are the required fields in a CST yaml settings file. The lists of frequencies specifies the frequency in each filename, so the lists must be in the same order (as must the feed_pol if it is a list):
telescope_name (str)
feed_name (str)
feed_version (str)
model_name (str)
model_version (str)
history (str)
frequencies (list(float))
cst text filenames (list(str)) – path relative to yaml file location
feed_pol (str) or (list(str))
Optional Fields
The following are optional fields:
ref_imp (float): beam model reference impedance
sim_beam_type (str): e.g. ‘E-farfield’
any other field that contains useful information that should be propagated with the beam. These will go into the extra_keywords attribute (note that if the field names are more than 8 characters they will be truncated to 8 if the beam is written to a beamfits file).
Example Settings File
An example settings yaml file for a HERA Vivaldi feed simulation is shown below. In this example, ‘software’, ‘layout’, and ‘port_num’ are extra fields that will be propagated with the beam. Note that quotes are used around version numbers that could be interpreted as floats and around filenames with spaces in them, but other strings do not require quotes. The history field shows a way to make a multi-line string (although this history is quite short, there could be many lines of history information there.):
telescope_name: HERA
feed_name: Vivaldi
feed_version: '1.0'
model_name: Mecha design - dish - cables - soil
model_version: '1.0'
software: CST 2016
history: |
beams simulated in Nov 2018 by NF
frequencies: [50e6, 51e6, 52e6, 53e6, 54e6, 55e6, 56e6, 57e6, 58e6, 59e6, 60e6, 61e6, 62e6, 63e6, 64e6, 65e6, 66e6, 67e6, 68e6, 69e6, 70e6, 71e6, 72e6, 73e6, 74e6, 75e6, 76e6, 77e6, 78e6, 79e6, 80e6, 81e6, 82e6, 83e6, 84e6, 85e6, 86e6, 87e6, 88e6, 89e6, 90e6, 91e6, 92e6, 93e6, 94e6, 95e6, 96e6, 97e6, 98e6, 99e6, 100e6, 101e6, 102e6, 103e6, 104e6, 105e6, 106e6, 107e6, 108e6, 109e6, 110e6, 111e6, 112e6, 113e6, 114e6, 115e6, 116e6, 117e6, 118e6, 119e6, 120e6, 121e6, 122e6, 123e6, 124e6, 125e6, 126e6, 127e6, 128e6, 129e6, 130e6, 131e6, 132e6, 133e6, 134e6, 135e6, 136e6, 137e6, 138e6, 139e6, 140e6, 141e6, 142e6, 143e6, 144e6, 145e6, 146e6, 147e6, 148e6, 149e6, 150e6, 151e6, 152e6, 153e6, 154e6, 155e6, 156e6, 157e6, 158e6, 159e6, 160e6, 161e6, 162e6, 163e6, 164e6, 165e6, 166e6, 167e6, 168e6, 169e6, 170e6, 171e6, 172e6, 173e6, 174e6, 175e6, 176e6, 177e6, 178e6, 179e6, 180e6, 181e6, 182e6, 183e6, 184e6, 185e6, 186e6, 187e6, 188e6, 189e6, 190e6, 191e6, 192e6, 193e6, 194e6, 195e6, 196e6, 197e6, 198e6, 199e6, 200e6, 201e6, 202e6, 203e6, 204e6, 205e6, 206e6, 207e6, 208e6, 209e6, 210e6, 211e6, 212e6, 213e6, 214e6, 215e6, 216e6, 217e6, 218e6, 219e6, 220e6, 221e6, 222e6, 223e6, 224e6, 225e6, 226e6, 227e6, 228e6, 229e6, 230e6, 231e6, 232e6, 233e6, 234e6, 235e6, 236e6, 237e6, 238e6, 239e6, 240e6, 241e6, 242e6, 243e6, 244e6, 245e6, 246e6, 247e6, 248e6, 249e6, 250e6]
filenames: ['farfield (f=50) [1].txt', 'farfield (f=51) [1].txt', 'farfield (f=52) [1].txt', 'farfield (f=53) [1].txt', 'farfield (f=54) [1].txt', 'farfield (f=55) [1].txt', 'farfield (f=56) [1].txt', 'farfield (f=57) [1].txt', 'farfield (f=58) [1].txt', 'farfield (f=59) [1].txt', 'farfield (f=60) [1].txt', 'farfield (f=61) [1].txt', 'farfield (f=62) [1].txt', 'farfield (f=63) [1].txt', 'farfield (f=64) [1].txt', 'farfield (f=65) [1].txt', 'farfield (f=66) [1].txt', 'farfield (f=67) [1].txt', 'farfield (f=68) [1].txt', 'farfield (f=69) [1].txt', 'farfield (f=70) [1].txt', 'farfield (f=71) [1].txt', 'farfield (f=72) [1].txt', 'farfield (f=73) [1].txt', 'farfield (f=74) [1].txt', 'farfield (f=75) [1].txt', 'farfield (f=76) [1].txt', 'farfield (f=77) [1].txt', 'farfield (f=78) [1].txt', 'farfield (f=79) [1].txt', 'farfield (f=80) [1].txt', 'farfield (f=81) [1].txt', 'farfield (f=82) [1].txt', 'farfield (f=83) [1].txt', 'farfield (f=84) [1].txt', 'farfield (f=85) [1].txt', 'farfield (f=86) [1].txt', 'farfield (f=87) [1].txt', 'farfield (f=88) [1].txt', 'farfield (f=89) [1].txt', 'farfield (f=90) [1].txt', 'farfield (f=91) [1].txt', 'farfield (f=92) [1].txt', 'farfield (f=93) [1].txt', 'farfield (f=94) [1].txt', 'farfield (f=95) [1].txt', 'farfield (f=96) [1].txt', 'farfield (f=97) [1].txt', 'farfield (f=98) [1].txt', 'farfield (f=99) [1].txt', 'farfield (f=100) [1].txt', 'farfield (f=101) [1].txt', 'farfield (f=102) [1].txt', 'farfield (f=103) [1].txt', 'farfield (f=104) [1].txt', 'farfield (f=105) [1].txt', 'farfield (f=106) [1].txt', 'farfield (f=107) [1].txt', 'farfield (f=108) [1].txt', 'farfield (f=109) [1].txt', 'farfield (f=110) [1].txt', 'farfield (f=111) [1].txt', 'farfield (f=112) [1].txt', 'farfield (f=113) [1].txt', 'farfield (f=114) [1].txt', 'farfield (f=115) [1].txt', 'farfield (f=116) [1].txt', 'farfield (f=117) [1].txt', 'farfield (f=118) [1].txt', 'farfield (f=119) [1].txt', 'farfield (f=120) [1].txt', 'farfield (f=121) [1].txt', 'farfield (f=122) [1].txt', 'farfield (f=123) [1].txt', 'farfield (f=124) [1].txt', 'farfield (f=125) [1].txt', 'farfield (f=126) [1].txt', 'farfield (f=127) [1].txt', 'farfield (f=128) [1].txt', 'farfield (f=129) [1].txt', 'farfield (f=130) [1].txt', 'farfield (f=131) [1].txt', 'farfield (f=132) [1].txt', 'farfield (f=133) [1].txt', 'farfield (f=134) [1].txt', 'farfield (f=135) [1].txt', 'farfield (f=136) [1].txt', 'farfield (f=137) [1].txt', 'farfield (f=138) [1].txt', 'farfield (f=139) [1].txt', 'farfield (f=140) [1].txt', 'farfield (f=141) [1].txt', 'farfield (f=142) [1].txt', 'farfield (f=143) [1].txt', 'farfield (f=144) [1].txt', 'farfield (f=145) [1].txt', 'farfield (f=146) [1].txt', 'farfield (f=147) [1].txt', 'farfield (f=148) [1].txt', 'farfield (f=149) [1].txt', 'farfield (f=150) [1].txt', 'farfield (f=151) [1].txt', 'farfield (f=152) [1].txt', 'farfield (f=153) [1].txt', 'farfield (f=154) [1].txt', 'farfield (f=155) [1].txt', 'farfield (f=156) [1].txt', 'farfield (f=157) [1].txt', 'farfield (f=158) [1].txt', 'farfield (f=159) [1].txt', 'farfield (f=160) [1].txt', 'farfield (f=161) [1].txt', 'farfield (f=162) [1].txt', 'farfield (f=163) [1].txt', 'farfield (f=164) [1].txt', 'farfield (f=165) [1].txt', 'farfield (f=166) [1].txt', 'farfield (f=167) [1].txt', 'farfield (f=168) [1].txt', 'farfield (f=169) [1].txt', 'farfield (f=170) [1].txt', 'farfield (f=171) [1].txt', 'farfield (f=172) [1].txt', 'farfield (f=173) [1].txt', 'farfield (f=174) [1].txt', 'farfield (f=175) [1].txt', 'farfield (f=176) [1].txt', 'farfield (f=177) [1].txt', 'farfield (f=178) [1].txt', 'farfield (f=179) [1].txt', 'farfield (f=180) [1].txt', 'farfield (f=181) [1].txt', 'farfield (f=182) [1].txt', 'farfield (f=183) [1].txt', 'farfield (f=184) [1].txt', 'farfield (f=185) [1].txt', 'farfield (f=186) [1].txt', 'farfield (f=187) [1].txt', 'farfield (f=188) [1].txt', 'farfield (f=189) [1].txt', 'farfield (f=190) [1].txt', 'farfield (f=191) [1].txt', 'farfield (f=192) [1].txt', 'farfield (f=193) [1].txt', 'farfield (f=194) [1].txt', 'farfield (f=195) [1].txt', 'farfield (f=196) [1].txt', 'farfield (f=197) [1].txt', 'farfield (f=198) [1].txt', 'farfield (f=199) [1].txt', 'farfield (f=200) [1].txt', 'farfield (f=201) [1].txt', 'farfield (f=202) [1].txt', 'farfield (f=203) [1].txt', 'farfield (f=204) [1].txt', 'farfield (f=205) [1].txt', 'farfield (f=206) [1].txt', 'farfield (f=207) [1].txt', 'farfield (f=208) [1].txt', 'farfield (f=209) [1].txt', 'farfield (f=210) [1].txt', 'farfield (f=211) [1].txt', 'farfield (f=212) [1].txt', 'farfield (f=213) [1].txt', 'farfield (f=214) [1].txt', 'farfield (f=215) [1].txt', 'farfield (f=216) [1].txt', 'farfield (f=217) [1].txt', 'farfield (f=218) [1].txt', 'farfield (f=219) [1].txt', 'farfield (f=220) [1].txt', 'farfield (f=221) [1].txt', 'farfield (f=222) [1].txt', 'farfield (f=223) [1].txt', 'farfield (f=224) [1].txt', 'farfield (f=225) [1].txt', 'farfield (f=226) [1].txt', 'farfield (f=227) [1].txt', 'farfield (f=228) [1].txt', 'farfield (f=229) [1].txt', 'farfield (f=230) [1].txt', 'farfield (f=231) [1].txt', 'farfield (f=232) [1].txt', 'farfield (f=233) [1].txt', 'farfield (f=234) [1].txt', 'farfield (f=235) [1].txt', 'farfield (f=236) [1].txt', 'farfield (f=237) [1].txt', 'farfield (f=238) [1].txt', 'farfield (f=239) [1].txt', 'farfield (f=240) [1].txt', 'farfield (f=241) [1].txt', 'farfield (f=242) [1].txt', 'farfield (f=243) [1].txt', 'farfield (f=244) [1].txt', 'farfield (f=245) [1].txt', 'farfield (f=246) [1].txt', 'farfield (f=247) [1].txt', 'farfield (f=248) [1].txt', 'farfield (f=249) [1].txt', 'farfield (f=250) [1].txt']
sim_beam_type: E-farfield
feed_pol: x
layout: 1 antenna
port_num: 1
ref_imp: 100
last updated: 2025-01-14