Full class list and documentation

Primary Beams

Simulations may be run using pyuvdata.UVBeam objects or pyuvdata.analytic_beam.AnalyticBeam objects. See the docs for those objects.

Antenna objects

class pyuvsim.Antenna(name: str, number: int, enu_position: ndarray[tuple[Any, ...], dtype[float]], beam_id: int)[source]

Describe a single interferometric element.

One of these defined for each antenna in the array.

Parameters
  • name (str) – Name of this antenna.

  • number (int) – Number of this antenna.

  • enu_position (array_like of float) – Position of this antenna in meters in the East, North, Up frame centered on the telescope location.

  • beam_id (int) – Index of the beam for this antenna from BeamList.

get_beam_jones(array, source_alt_az, frequency, reuse_spline=True, interpolation_function=None, freq_interp_kind=None, beam_interp_check=True)[source]

Calculate the jones matrix for this antenna in the direction of sources.

A Nfeeds x 2 x Nsources array of Efield vectors in Az/Alt.

Parameters
  • array (Telescope object) – Provides the beam list.

  • source_alt_az (array_like) – Positions to evaluate in alt/az, where source_alt_az[0] gives list of alts soruce_alt_az[1] gives list of corresponding az

  • frequency (float or Quantity) – Frequency. Assumed to be Hz if float.

  • reuse_spline (bool) – Option to keep and reuse interpolation splines in pyuvdata.UVBeam.

  • interpolation_function (str) – Set the angular interpolation function on the pyuvdata.UVBeam. See pyuvdata.UVBeam.interp() for options.

  • freq_interp_kind (str) – Interpolation method for frequencies. See pyuvdata.UVBeam.interp() for options.

  • beam_interp_check (bool) – Option to check whether the beam covers all the skymodel az/za values to ensure that they are covered by the intrinsic data array. Checking them can be quite computationally expensive. Defaults to True to run the check, can be set to False if the beam coverage is known.

Returns

jones_matrix (array_like of float) – Jones matricies for each source location, shape (Nfeeds,2, Ncomponents). The first axis is feed, the second axis is vector component on the sky in az/za.

Baseline objects

class pyuvsim.Baseline(antenna1, antenna2)[source]

Defines a single pair of Antenna objects with useful methods.

Defines the means for calculating their uvw coordinates, and comparing them.

Parameters

antenna1, antenna2 (pyuvsim.Antenna) – The two antennas that make up this baseline.

Attributes
  • antenna1, antenna2 (pyuvsim.Antenna) – The two antennas that make up this baseline.

  • enu (array_like of float) – The vector pointing from antenna1 to antenna2 in ENU coordinates (relative to the telescope location.)

  • uvw (array_like of float) – The uvw coordinate of this baseline. Identical to enu since we work in the local alt/az frame.

Telescopes

Shared properties of all antennas, such as array center location and list of primary beam objects.

class pyuvsim.Telescope(telescope_name: str, telescope_location: EarthLocation, beam_list: BeamList)[source]

Container for data common to all antennas in the array.

Defines the location and name of the observing site, and holds the list of beam objects used by the array.

Parameters
class pyuvsim.BeamList(beam_list: [pyuvdata.uvbeam.uvbeam.UVBeam | pyuvdata.analytic_beam.AnalyticBeam], *, beam_type: Optional[Literal['efield', 'power']] = 'efield', spline_interp_opts: dict[str, int] | None = None, freq_interp_kind: str = 'cubic', peak_normalize: dataclasses.InitVar[bool] = True)[source]

A container for the set of beam models and related parameters.

When used with MPI, this can be initialized simultaneously on all processes such that sky_in is provided only on the root process. The data can then be shared across all processes by running the share method.

This behaves just as a normal Python list in terms of indexing but it does not allow changes to the list.

Attributes
  • beam_list (list of pyuvdata.BeamInterface) – The list of BeamInterface objects.

  • spline_interp_opts (dict) – Degrees of the bivariate spline for the angular interpolation (passed directly to numpy.RectBivariateSpline, e.g. use {‘kx’ : 2, ‘ky’ : 2} for 2nd degree splines in each direction.) Default is 3rd order in each direction.

  • freq_interp_kind (str) – The type of frequency interpolation, anything accepted by scipy.interpolate.interp1d. Default is “cubic”.

Parameters
  • beam_list (list pyuvdata.UVBeam or pyuvdata.AnalyticBeam) – A list of pyuvdata UVBeam or AnalyticBeam objects. These will be converted into BeamInterface objects during initialization.

  • beam_type (str) – Desired beam type, either “efield” or “power”, defaults to “efield”.

  • spline_interp_opts (dict) – Provide options to numpy.RectBivariateSpline. This includes spline order parameters kx and ky, and smoothing parameter s. Only applies UVBeams and for az_za_simple interpolation.

  • freq_interp_kind (str) – Interpolation method to use for frequencies. See scipy.interpolate.interp1d for details. Defaults to “cubic”.

  • peak_normalize (bool) – Option to peak normalize UVBeams in the beam list. This is required to be True for the pyuvsim simulator. Defaults to True.

check_all_azza_beams_full_sky()[source]

Check if all az_za UVBeams cover the full sky.

Used to decide whether to turn off checking that the beam covers the interpolation location.

property data_normalization

Return the data_normalization of all beams in list.

share(root=0)[source]

Share across MPI processes (requires mpi4py to use).

All attributes are put in shared memory.

Parameters

root (int) – Root rank on COMM_WORLD, from which data will be broadcast.

Sources

See documentation for pyradiosky.SkyModel.

Simulation setup

Simulations are run from .yaml files with specified keywords. See Running a simulation.

Functions for setting up simulations based on input config files.

This provides a configuration interface of human-readable yaml and csv files to specify all the required simulation parameters. This configuration interface is used by other simulators as well.

This module contains methods to create configuration files from pyuvdata.UVData objects and empty pyuvdata.UVData objects from configuration files.

class pyuvsim.simsetup.SkyModelData(sky_in=None, filename=None)[source]

Carries immutable SkyModel data in simple ndarrays.

This is to facilitate sharing SkyModel objects in MPI, without excessive copying and memory bloat.

When used with MPI, this can be initialized simultaneously on all processes such that sky_in is provided only on the root process. The data can then be shared across all processes by running the share method.

Parameters
  • sky_in (pyradiosky.SkyModel) – A valid SkyModel object.

  • filename (str or list of str, optional) – The filename (or other string identifier) of the input catalog. This overrides the filename set on the sky_in object (if it has one). If not set, this defaults to the filename set on the sky_in object (if it has one)

get_skymodel(inds=None)[source]

Initialize pyradiosky.SkyModel from current settings.

Parameters

inds (range or index array) – Indices to select along the component axis.

share(root=0)[source]

Share across MPI processes (requires mpi4py to use).

All attributes are put in shared memory.

Parameters

root (int) – Root rank on COMM_WORLD, from which data will be broadcast.

subselect(inds)[source]

Subselect, returning a new SkyModelData object.

Parameters

inds (range or index array) – Indices to select along the component axis.

Returns

SkyModelData – A new SkyModelData with Ncomp axes downselected.

Notes

If inds is a range object, this method will avoid copying data in numpy arrays, such that the returned SkyModelData object carries views into the current object’s arrays.

pyuvsim.simsetup.create_mock_catalog(time, arrangement='zenith', array_location=None, Nsrcs=None, alt=None, save=False, min_alt=None, rseed=None, return_data=False, diffuse_model=None, diffuse_params=None, map_nside=None)[source]

Create a mock catalog.

SkyModels are defined in an AltAz frame at the given time, then returned in ICRS ra/dec coordinates. They always have a ‘flat’ spectral type.

Parameters
  • time (float or astropy.time.Time object) – Julian date, if a float it is interpreted as a JD.

  • arrangement (str) – Point source pattern (default = 1 source at zenith). Accepted arrangements:

    • triangle: Three point sources forming a triangle around the zenith

    • cross: An asymmetric cross

    • zenith: Some number of sources placed at the zenith.

    • off-zenith: A single source off zenith

    • long-line: Horizon to horizon line of point sources

    • hera_text: Spell out HERA around the zenith

    • random: Randomly distributed point sources near zenith

    • diffuse: Choice of solved models in the analytic_diffuse module.

  • Nsrcs (int) – Number of sources to make (used for zenith, off-zenith, long-line, and random arrangements).

  • array_location (astropy.coordinates.EarthLocation) – Location of observer. Source positions are defined with respect to a particular zenith. Can also be a lunarsky.MoonLocation object to make catalogs for observers on the Moon. Default = HERA site

  • alt (float) – For off-zenith and triangle arrangements, altitude to place sources. In degrees.

  • min_alt (float) – For random and long-line arrangements, minimum altitude at which to place sources. In degrees.

  • save (bool) – Save mock catalog as npz file.

  • rseed (int) – If using the random configuration, pass in a RandomState seed.

  • return_data (bool) – If True, return a SkyModelData object instead of pyradiosky.SkyModel.

  • diffuse_model (str) – If arrangement is ‘diffuse’, name of the diffuse model to generate. See documentation in analytic_diffuse.models for details.

  • diffuse_params (dict) – If arrangement is ‘diffuse’, extra parameters accepted by the chosen model. See documentation in analytic_diffuse.models for details.

  • map_nside (int) – If arrangement is ‘diffuse’, nside parameter for the healpix map to make.

Returns

  • pyradiosky.SkyModel or SkyModelData – The catalog, as either a SkyModel or a SkyModelData (if return_data is True)

  • dict – A dictionary of keywords used to define the catalog.

Notes

Generating diffuse models requires the analytic_diffuse and astropy_healpix modules.

pyuvsim.simsetup.initialize_catalog_from_params(obs_params, input_uv=None, filetype=None, return_catname=False)[source]

Make catalog from parameter file specifications.

Default behavior is to do coarse horizon cuts.

Parameters
  • obs_params (str or dict) – Either an obsparam file name or a dictionary of parameters.

  • input_uv (pyuvdata.UVData) – Used to set location for mock catalogs and for horizon cuts.

  • filetype (str) – One of [‘skyh5’, ‘gleam’, ‘vot’, ‘text’, ‘fhd’] or None. If None, the code attempts to guess what the file type is.

  • return_catname (bool) – Return the catalog name. Default is False.

Returns

  • sky (pyradiosky.SkyModel) – Source catalog filled with data.

  • source_list_name (str) – Catalog identifier for metadata. Only returned if return_catname is True.

pyuvsim.simsetup.initialize_uvdata_from_keywords(output_yaml_filename=None, antenna_layout_filepath=None, output_layout_filename=None, array_layout=None, telescope_location=None, telescope_name=None, feed_array=None, feed_angle=None, mount_type=None, Nfreqs=None, start_freq=None, bandwidth=None, freq_array=None, channel_width=None, Ntimes=None, integration_time=None, start_time=None, time_array=None, bls=None, antenna_nums=None, antenna_names=None, polarization_array=None, no_autos=False, redundant_threshold=None, conjugation_convention=None, blt_order=None, write_files=True, path_out=None, complete=False, check_kw: dict | None = None, **kwargs)[source]

Initialize a pyuvdata.UVData object from keyword arguments.

Optionally, write out the configuration to YAML and CSV files such that initialize_uvdata_from_params() will produce the same pyuvdata.UVData object.

Parameters
  • output_yaml_filename (str (optional)) – Specify filename for yaml file to write out. Defaults to obsparam.yaml

  • antenna_layout_filepath (str (optional)) – Path to csv file of antenna positions (see documentation for details).

  • output_layout_filename (str (optional)) – File name for antenna layout if writing files out. If unspecified, defaults to the basename of the antenna_layout_filepath. If that is not given, it defaults to “antenna_layout.csv” Will not overwrite existing files.

  • array_layout (dictionary (required if antenna_layout_filepath not given)) – keys are integer antenna numbers, values are len-3 antenna positions in ENU coordinates [meters].

  • antenna_names (list of str (optional)) – If unset, antenna names are assigned as “%s” % antnum.

  • telescope_location (array of float, shape (3,)) – Telescope location on Earth in LatLonAlt coordinates [deg, deg, meters]

  • telescope_name (str) – Name of telescope

  • feed_array (array-like of str or None) – List of feeds for each antenna in the telescope, must be one of “x”, “y”, “l”, “r”. Shape (Nants, Nfeeds), dtype str. Can also be shape (Nfeeds,), in which case the same values for feed_array are used for all antennas in the object.

  • feed_angle (array-like of float or None) – Orientation of the feed with respect to zenith (or with respect to north if pointed at zenith). Units is in rads, vertical polarization is nominally 0, and horizontal polarization is nominally pi / 2. Shape (Nants, Nfeeds), dtype float. Can also be shape (Nfeeds,), in which case the same values for feed_angle are used for all antennas in the object.

  • mount_type (str or array-like of 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. Shape (Nants,), dtype str. Can also provide a single string, in which case the same mount_type is used for all antennas in the object.

  • Nfreqs (int) – Number of frequency channels

  • start_freq (float) – Starting frequency [Hz]

  • bandwidth (float) – Total frequency bandwidth of spectral window [Hz]

  • channel_width (float) – Frequency channel spacing [Hz]

  • freq_array (ndarray) – frequency array [Hz], if this is specified, it supersedes Nfreqs, start_freq, bandwidth

  • Ntimes (int) – Number of integration bins

  • integration_time (float) – Width of time bins [seconds]

  • start_time (float) – Time of the first integration bin [Julian Date]

  • time_array (ndarray) – time array [Julian Date]. If this is specified it supersedes values of Ntimes, start_time and integration_time

  • bls (list) – List of antenna-pair tuples for baseline selection

  • redundant_threshold (float) – Redundant baseline selection tolerance for selection [meters]

  • antenna_nums (list) – List of antenna numbers to keep in array

  • polarization_array (list) – List of polarization strings (or ints) to insert into object

  • no_autos (bool) – If True, eliminate all auto correlations

  • conjugation_convention (str, optional) – A convention for the directions of the baselines, options are: ‘ant1<ant2’, ‘ant2<ant1’, ‘u<0’, ‘u>0’, ‘v<0’, ‘v>0’

  • blt_order (str or array-like of str, optional) – Specifies the ordering along the blt axis. Defaults to [‘time’, ‘baseline’]. If set to anything else, it is passed to pyuvdata.UVData.reorder_blts(). A single string specifies the order, a list of two strings specifies the [order, minor_order]. See the docs on the UVData method for more information.

  • write_files (bool) – If True, write out the parameter information to yaml files.

  • path_out (str (optional)) – Path in which to place generated configuration files, if write_files is True. Defaults to current directory.

  • complete (bool (optional)) – Whether to fill out the pyuvdata.UVData object with its requisite data arrays, and check if it’s all consistent.

  • check_kw (dict (optional)) – A dictionary of keyword arguments to pass to the uvdata.UVData.check() method after object creation. Caution: turning off critical checks can result in a UVData object that cannot be written to a file.

  • kwargs (dictionary) – Any additional valid pyuvdata.UVData attribute to assign to object.

Returns

pyuvdata.UVData – Initialized based on keywords, with a zeroed data_array, no flags and nsample_array of all ones.

pyuvsim.simsetup.initialize_uvdata_from_params(obs_params, return_beams=False, reorder_blt_kw=None, check_kw=None)[source]

Construct a pyuvdata.UVData object from parameters in a valid yaml file.

Sufficient information must be provided by the parameters to define time and frequency arrays and verify the channel widths and time steps. This will error if insufficient or incompatible parameters are defined.

The parameter dictionary may contain any valid pyuvdata.UVData attributes as well.

If the polarization array is not specified, it defaults to (XX, XY, YX, YY).

Parameters
  • obs_params (dict or str) – Either an obs_param file name or a dictionary of parameters read in. Additional pyuvdata.UVData parameters may be passed in through here.

  • return_beams (bool) – Option to return the beam_list and beam_dict. Default is False.

  • reorder_blt_kw (dict (optional)) – Deprecated. Use obs_param[‘ordering’][‘blt_order’] instead. Keyword arguments to send to the uvdata.reorder_blts method at the end of the setup process. Typical parameters include “order” and “minor_order”. Default values are order='time' and minor_order='baseline'.

  • check_kw (dict (optional)) – A dictionary of keyword arguments to pass to the uvdata.UVData.check() method after object creation. Caution: turning off critical checks can result in a UVData object that cannot be written to a file.

Returns

  • uv_obj (pyuvdata.UVData) – Initialized UVData object.

  • beam_list (pyuvsim.BeamList) – List of beam specifiers as strings, only returned if return_beams is True.

  • beam_dict (dict) – Map of antenna numbers to index in beam_list, only returned if return_beams is True.

pyuvsim.simsetup.uvdata_to_config_file(uvdata_in, param_filename=None, telescope_config_name='', layout_csv_name='', catalog='mock', path_out='.')[source]

Extract simulation configuration settings from a UVData object.

When used with uvdata_to_telescope_config(), this will produce all the necessary configuration yaml and csv file to make an “empty” pyuvdata.UVData object comparable to the argument uvdata_in. The generated file will match uvdata_in in frequency, time, antenna positions, and uvw coordinates. Note that it may be different in the baseline conjugation convention and the ordering of the baseline-time axis, but those can be modified using the pyuvdata.UVData.conjugate_bls() and pyuvdata.UVData.reorder_blt() methods.

Parameters
  • uvdata_in (pyuvdata.UVData) – UVData object for which to make config files.

  • param_filename (str) – output param file name, defaults to obsparam_#.yaml.

  • telescope_config_name (str) – Name of telescope configuration yaml file. Defaults to blank string.

  • layout_csv_name (str) – Name of antenna layout csv file. Defaults to blank string.

  • catalog (str) – Path to catalog file. Defaults to ‘mock’.

  • path_out (str) – Where to put config files.

pyuvsim.simsetup.uvdata_to_telescope_config(uvdata_in, beam_filepath, layout_csv_name=None, telescope_config_name=None, return_names=False, path_out='.')[source]

Make telescope parameter files from a pyuvdata.UVData object.

Makes both a telescope_config file and a layout_csv file:

  • telescope_config: YAML file with telescope_location and telescope_name

    The beam list is spoofed, since that information cannot be found in a UVData object.

  • layout_csv: tab separated value file giving ENU antenna positions.

    Beam ID is spoofed as well.

See https://pyuvsim.readthedocs.io/en/latest/parameter_files.html for details.

Parameters
  • uvdata_in (pyuvdata.UVData) – UVData object for which to make config files.

  • path_out (str) – Target directory for the config files.

  • beam_filepath (str) – Path to a beamfits file.

  • layout_csv_name (str) – The name for the antenna positions. Default <telescope_name>_layout.csv, where <telescope_name> is uvdata_in.telescope_name

  • telescope_config_name (str) – The name for the telescope config file Default telescope_config_<telescope_name>.yaml

  • return_names (bool) – Return the file names. Default is False. Used in tests.

Returns

tuple – if return_names, returns (path, telescope_config_name, layout_csv_name)

UV Simulation functions

Methods for running simulations.

The primary machinery of pyuvsim.

The UVTask and UVEngine classes and the functions that actually run the simulation.

class pyuvsim.uvsim.UVEngine(task=None, update_positions=True, update_beams=True, reuse_spline=True)[source]

The object where the visibility calculations actually take place.

Parameters
  • task (UVTask) – The task to do the calculations for.

  • reuse_spline (bool) – Option to reuse the spline in the beam interpolation to save time.

  • update_positions (bool) – Flag indicating that positions need to be updated (when time or source changes).

  • update beams (bool) – Flag indicating that beams need to be updated (when time, sources, frequency, or beam pair changes).

Attributes
  • reuse_spline (bool) – Option to reuse the spline in the beam interpolation to save time.

  • update_positions (bool) – Flag indicating that positions need to be updated (when time or source changes).

  • update_local_coherency (bool) – Flag indicating that local coherencies need to be updated (when time or source changes).

  • update beams (bool) – Flag indicating that beams need to be updated (when time, sources, frequency, or beam pair changes).

  • sources (pyuvsim.simsetup.SkyModelData) – The sources to include in the visibility.

  • current_time (astropy.time.Time) – Time for the current calculation.

  • current_freq (astropy.units.Quantity) – Frequency for the current calculation, units compatible with Hz.

  • current_beam_pair (2-tuple of ints) – Tuple containing the beam ids for the current calculation.

  • beam1_jones (array_like of float) – Jones matricies for the first beam for each source location, shape (2,2, Ncomponents). The first axis is feed, the second axis is vector component on the sky in az/za.

  • beam2_jones (array_like of float) – Jones matricies for the second beam for each source location, shape (2,2, Ncomponents). The first axis is feed, the second axis is vector component on the sky in az/za.

  • local_coherency (array_like of float) – Source local coherencies in alt/az basis, shape (2, 2, Nfreqs, Ncomponents).

  • apparent_coherency (array_like of float) – Source apparent coherencies (including the beam response), shape (2, 2, Nfreqs, Ncomponents).

  • task (UVTask) – The task currently being calculated.

apply_beam(beam_interp_check=True)[source]

Set apparent coherency from jones matrices and source coherency.

beam_interp_checkbool

Option to enable checking that the source positions are within the area covered by the beam. If all the az/za beams cover the full sky horizon to horizon this checking is turned off by default in run_uvdata_uvsim, otherwise it is turned on. Setting this to False can speed up simulations but if sources are simulated outside the beam area the response will be incorrect. This keyword only applies to beams that are regularly gridded in azimuth and zenith angle.

make_visibility()[source]

Calculate visibility contribution from a set of source components.

set_task(task)[source]

Set the task for the Engine.

Parameters

task (UVTask) – The task to do the calculations for.

class pyuvsim.uvsim.UVTask(sources, time, freq, baseline, telescope, freq_i=0)[source]

An object to hold all the information to calculate a visibility for a set of sources.

Parameters
  • sources (pyradiosky.SkyModel) – The sources to include in the visibility.

  • time (astropy.time.Time object or float) – Time at which to calculate the visibility, either an astropy Time object or a float in Julian Day (will be converted to a Time object using format=’jd’).

  • freq (astropy.units.Quantity or float) – Frequency at which to calculate the visibility, either an astropy Quantity with units compatible with Hz or a float in Hz (will be converted to an astropy Quantity with units of Hz)

  • baseline (pyuvsim.Baseline) – The baseline to calculate the visibility for.

  • telescope (pyuvsim.Telescope) – The telescope object this visibility belongs to, which carries the telescope location and beam information.

  • freq_i (int) – Frequency index for this visibility’s source coherency.

Attributes
  • sources (pyradiosky.SkyModel) – The sources to include in the visibility.

  • time (astropy.time.Time) – Time at which to calculate the visibility.

  • freq (astropy.units.Quantity) – Frequency at which to calculate the visibility.

  • baseline (pyuvsim.Baseline) – The baseline to calculate the visibility for.

  • telescope (pyuvsim.Telescope) – The telescope object this visibility belongs to.

  • freq_i (int) – Frequency index for this visibility’s source coherency. Always 0 for flat spectrum sources, otherwise should be the index of freq in the full frequency array for this simulation.

  • uvdata_index (tuple of int) – Indexes for this visibility location in the uvdata data_array, length 3 giving (blt_index, 0, frequency_index), where the zero index is for the old spw axis.

  • visibility_vector (ndarray of complex) – The calculated visibility, shape (Npols,) ordered as [xx, yy, xy, yx].

pyuvsim.uvsim.run_uvdata_uvsim(input_uv: UVData, beam_list: BeamList, beam_dict: dict, catalog: pyradiosky.skymodel.SkyModel | pyuvsim.simsetup.SkyModelData, beam_interp_check: bool | None = None, quiet: bool = False, block_nonroot_stdout: bool = True, Nsky_parts: int | None = None) UVData[source]

Run uvsim from UVData object.

Parameters
  • input_uv (pyuvdata.UVData instance) – Provides baseline/time/frequency information.

  • beam_list (pyuvsim.BeamList) – BeamList carrying beam models.

  • beam_dict (dictionary, optional) – {antenna_name : beam_id}, where beam_id is an index in the beam_list. This is used to assign beams to antennas. Default: All antennas get the 0th beam in the beam_list.

  • catalog (pyuvsim.simsetup.SkyModelData or pyradiosky.SkyModel) – Source catalog, either a pyradiosky SkyModel or a SkyModelData object.

  • beam_interp_check (bool) – Option to enable checking that the source positions are within the area covered by the beam. If the beam covers the full sky horizon to horizon this checking is turned off by default, otherwise it is turned on. Setting this to False can speed up simulations but if sources are simulated outside the beam area the response will be incorrect. This keyword only applies to beams that are regularly gridded in azimuth and zenith angle.

  • quiet (bool) – Do not print anything.

  • block_nonroot_stdout (bool) – Redirect stdout on nonzero ranks to /dev/null, for cleaner output.

  • Nsky_parts (int, optional) – Number of parts to chunk the source list into to reduce memory. The default is to use the minimum number to fit within the memory of the processing unit. This parameter is included mostly to allow for testing, but can be used by a knowledgeable user if needed. If set too low, a ValueError will be raised.

Returns

pyuvdata.UVData – instance containing simulated visibilities.

pyuvsim.uvsim.run_uvsim(params, return_uv=False, beam_interp_check=None, quiet=False, block_nonroot_stdout=True)[source]

Run a simulation off of an obsparam yaml file.

Parameters
  • params (str) – Path to a parameter yaml file or a parameter dict.

  • return_uv (bool) – If true, do not write results to file and return uv_out. (Default False)

  • beam_interp_check (bool) – Option to enable checking that the source positions are within the area covered by the beam. If the beam covers the full sky horizon to horizon this checking is turned off by default, otherwise it is turned on. Setting this to False can speed up simulations but if sources are simulated outside the beam area the response will be incorrect. This keyword only applies to beams that are regularly gridded in azimuth and zenith angle.

  • quiet (bool) – If True, do not print anything to stdout. (Default False)

  • block_nonroot_stdout (bool) – Redirect stdout on nonzero ranks to /dev/null, for cleaner output.

Returns

uv_out (pyuvdata.UVData) – Finished simulation results. Returned only if return_uv is True.

pyuvsim.uvsim.uvdata_to_task_iter(task_ids, input_uv, catalog, beam_list, beam_dict, Nsky_parts=1)[source]

Generate UVTask objects.

Parameters
  • task_ids (range) – Task indices in the full flattened meshgrid of parameters.

  • input_uv (pyuvdata.UVData) – UVData object to be filled with data.

  • catalog (pyuvsim.simsetup.SkyModelData) – Source components.

  • beam_list (pyuvsim.BeamList) – BeamList carrying beam model objects.

  • beam_dict (dict) – Map of antenna numbers to index in beam_list.

  • Nsky_parts (int) – Number of chunks to break the skymodel into. Default is 1.

Yields

Iterable of UVTask objects.

MPI Tools

Tools to initialize MPI and share data among processors.

Profiling

Developers will want to use line profiling tools to track down bottlenecks and find room for performance improvements. These tools better allow profilers to operate within the MPI environment.

Use the line profiler when requested.

pyuvsim.profiling.LineProfiler()[source]

Mock to fix imports.

pyuvsim.profiling.get_profiler()[source]

Get the profiler.

pyuvsim.profiling.set_profiler(func_list=['get_beam_jones', 'initialize_uvdata_from_params', 'apply_beam', 'make_visibility', 'update_positions', 'coherency_calc', 'uvdata_to_task_iter', 'run_uvdata_uvsim', 'run_uvsim'], rank=0, outfile_prefix='time_profile.out', dump_raw=False)[source]

Apply a line profiler to the listed functions, wherever they appear in pyuvsim.

Places a LineProfiler object in the module namespace, and registers its dumping/printing functions to run at the end. When the Python environment closes, the profiler functions print_stats (and dump_stats, if dump_raw is True) will execute, saving profiler data to file.

Parameters
  • func_list (list) – List of function names (strings) to profile. Defaults to profiling.default_profile_funcs.

  • rank (int, optional) – Which rank process should write out to file? (only one rank at a time will).

  • outfile_prefix (str) –

    Filename prefix for printing profiling results.

    Human-readable line by line profiling goes to <outfile_prefix>.out LineStats data goes to <outfile_prefix>.lprof (if dump_raw) Axis sizes go to <outfile_prefix>_axes.npz

  • dump_raw (bool) – Write out a pickled LineStats object to <outfile_name>.lprof (Default False)