karabo.simulation

class SkyModel(sources: ~numpy.ndarray[~typing.Any, ~numpy.dtype[~numpy.float64]] | ~numpy.ndarray[~typing.Any, ~numpy.dtype[~numpy.object_]] | ~xarray.core.dataarray.DataArray | None = None, wcs: ~astropy.wcs.wcs.WCS | None = None, precision: ~typing.Type[~numpy.float64] = <class 'numpy.float64'>, h5_file_connection: ~h5py._hl.files.File | None = None)

Class containing all information of the to be observed Sky.

SkyModel.sources is a xarray.DataArray

( https://docs.xarray.dev/en/stable/generated/xarray.DataArray.html ).

np.ndarray are also supported as input type for SkyModel.sources,

however, the values in SkyModel.sources are converted to xarray.DataArray.

SkyModel.compute method is used to load the data into memory as a numpy array. It should be called after all the filtering and other operations are completed and to avoid doing the same calculation multiple them when e.g. on a cluster.

Variables:
  • sources

    List of all point sources in the sky as xarray.DataArray. The source_ids reside in SkyModel.source_ids if provided through xarray.sources.coords with an arbitrary string key as index or np.ndarray as idx SOURCES_COLS. A single point source is described using the following col-order:

    • [0] right ascension (deg)

    • [1] declination (deg)

    • [2] stokes I Flux (Jy)

    • [3] stokes Q Flux (Jy): defaults to 0

    • [4] stokes U Flux (Jy): defaults to 0

    • [5] stokes V Flux (Jy): defaults to 0

    • [6] reference_frequency (Hz): defaults to 0

    • [7] spectral index (N/A): defaults to 0

    • [8] rotation measure (rad / m^2): defaults to 0

    • [9] major axis FWHM (arcsec): defaults to 0

    • [10] minor axis FWHM (arcsec): defaults to 0

    • [11] position angle (deg): defaults to 0

    • [12] true redshift: defaults to 0

    • [13] observed redshift: defaults to 0

    • [14] object-id: just for np.ndarray

      it is removed in the xr.DataArray and exists then in xr.DataArray.coords as index.

  • wcs – World Coordinate System (WCS) object representing the coordinate transformation between pixel coordinates and celestial coordinates (e.g., right ascension and declination).

  • precision – The precision of numerical values used in the SkyModel. Has to be of type np.float_.

  • h5_file_connection – An open connection to an HDF5 (h5) file that can be used to store or retrieve data related to the SkyModel.

__init__(sources: ~numpy.ndarray[~typing.Any, ~numpy.dtype[~numpy.float64]] | ~numpy.ndarray[~typing.Any, ~numpy.dtype[~numpy.object_]] | ~xarray.core.dataarray.DataArray | None = None, wcs: ~astropy.wcs.wcs.WCS | None = None, precision: ~typing.Type[~numpy.float64] = <class 'numpy.float64'>, h5_file_connection: ~h5py._hl.files.File | None = None) None

Initialize a SkyModel object.

Parameters

sources{xarray.DataArray, np.ndarray}, optional

List of all point sources in the sky. It can be provided as an xarray.DataArray or np.ndarray. If provided as an np.ndarray, the values are converted to xarray.DataArray.

wcsWCS, optional

World Coordinate System (WCS) object representing the coordinate transformation between pixel coordinates and celestial coordinates.

precisionnp.dtype, optional

The precision of numerical values used in the SkyModel. It should be a NumPy data type (e.g., np.float64).

h5_file_connectionh5py.File, optional

An open connection to an HDF5 (h5) file that can be used to store or retrieve data related to the SkyModel.

add_point_sources(sources: ndarray[Any, dtype[float64]] | ndarray[Any, dtype[object_]] | DataArray) None

Add new point sources to the sky model.

Parameters:

sourcesnp.ndarray with shape (number of sources, 1 + SOURCES_COLS),

where you can place the “source_id” at index SOURCES_COLS. OR an xarray.DataArray with shape (number of sources, SOURCES_COLS), where you can place the “source_id” at xarray.DataArray.coord or use SkyModel.source_ids later.

The column indices correspond to:

  • [0] right ascension (deg)

  • [1] declination (deg)

  • [2] stokes I Flux (Jy)

  • [3] stokes Q Flux (Jy): defaults to 0

  • [4] stokes U Flux (Jy): defaults to 0

  • [5] stokes V Flux (Jy): defaults to 0

  • [6] reference_frequency (Hz): defaults to 0

  • [7] spectral index (N/A): defaults to 0

  • [8] rotation measure (rad / m^2): defaults to 0

  • [9] major axis FWHM (arcsec): defaults to 0

  • [10] minor axis FWHM (arcsec): defaults to 0

  • [11] position angle (deg): defaults to 0

  • [12] true redshift: defaults to 0

  • [13] observed redshift: defaults to 0

  • [14] source id (object): is in SkyModel.source_ids if provided

close() None

Closes the connection to the HDF5 file.

This method closes the connection to the HDF5 file if it is open and sets the h5_file_connection attribute to None.

compute() None

Loads the lazy data into a numpy array, wrapped in a xarray.DataArray.

This method loads the lazy data from the sources into a numpy array, which is then wrapped in a xarray.DataArray object. It performs the computation necessary to obtain the actual data and stores it in the _sources attribute. After the computation is complete, it calls the close method to close the connection to the HDF5 file.

Returns: None

convert_to_backend(backend: Literal[SimulatorBackend.OSKAR] = SimulatorBackend.OSKAR, desired_frequencies_hz: Literal[None] = None, channel_bandwidth_hz: float | None = None, verbose: bool = False) SkyModel
convert_to_backend(backend: Literal[SimulatorBackend.RASCIL], desired_frequencies_hz: ndarray[Any, dtype[float64]], channel_bandwidth_hz: float | None = None, verbose: bool = False) List[SkyComponent]

Convert an existing SkyModel instance into a format acceptable by a desired backend. backend: Determines how to return the SkyModel source catalog.

OSKAR: return the current SkyModel instance, since methods in Karabo support OSKAR-formatted source np.array values. RASCIL: convert the current source array into a list of RASCIL SkyComponent instances.

desired_frequencies_hz: List of frequencies corresponding to start of desired frequency channels. This field is required to convert sources into RASCIL SkyComponents.

The array contains starting frequencies for the desired channels. E.g. [100e6, 110e6] corresponds to 2 frequency channels, which start at 100 MHz and 110 MHz, both with a bandwidth of 10 MHz.

channel_bandwidth_hz: Used if desired_frequencies_hz has only one element.

Otherwise, bandwidth is determined as the delta between the first two entries in desired_frequencies_hz.

verbose: Determines whether to display additional print statements.

explore_sky(phase_center: ~typing.List[int] | ~typing.List[float], stokes: ~typing.Literal['Stokes I', 'Stokes Q', 'Stokes U', 'Stokes V'] = 'Stokes I', idx_to_plot: ~numpy.ndarray[~typing.Any, ~numpy.dtype[~numpy.int64]] | None = None, xlim: ~typing.Tuple[int | float, int | float] | None = None, ylim: ~typing.Tuple[int | float, int | float] | None = None, figsize: ~typing.Tuple[int | float, int | float] | None = None, title: str | None = None, xlabel: str | None = None, ylabel: str | None = None, cfun: ~typing.Callable[[...], bool | int | ~numpy.integer | float | ~numpy.floating | ~numpy.ndarray[~typing.Any, ~numpy.dtype[~numpy.int64 | ~numpy.float64]]] | None = <ufunc 'log10'>, cmap: str | None = 'plasma', cbar_label: str | None = None, with_labels: bool = False, wcs: ~astropy.wcs.wcs.WCS | None = None, wcs_enabled: bool = True, filename: str | None = None, block: bool = False, **kwargs: ~typing.Any) None

A scatter plot of the SkyModel (self) where the sources are projected according to the phase_center

Parameters:
  • phase_center – [RA,DEC]

  • stokesSkyModel stoke flux

  • idx_to_plot – If you want to plot only a subset of the sources, set the indices here.

  • xlim – RA-limit of plot

  • ylim – DEC-limit of plot

  • figsize – figsize as tuple

  • title – plot title

  • xlabel – xlabel

  • ylabel – ylabel

  • cfun – flux scale transformation function for scatter-coloring

  • cmap – matplotlib color map

  • cbar_label – color bar label

  • with_labels – Plots object ID’s if set?

  • wcs – If you want to use a custom astropy.wcs, ignores phase_center if set

  • wcs_enabled – Use wcs transformation?

  • filename – Set to path/fname to save figure (set extension to fname to overwrite .png default)

  • block – Whether or not plotting should block the rest of the program

  • kwargs – matplotlib kwargs for scatter & Collections, e.g. customize s, vmin or vmax

filter_by_column(col_idx: int, min_val: int | float, max_val: int | float) SkyModel

Filters the sky based on a specific column index

Parameters:
  • col_idx – Column index to filter by

  • min_val – Minimum value for the column

  • max_val – Maximum value for the column

Return sky:

Filtered copy of the sky

filter_by_radius(inner_radius_deg: int | float, outer_radius_deg: int | float, ra0_deg: int | float, dec0_deg: int | float, indices: Literal[False] = False) SkyModel
filter_by_radius(inner_radius_deg: int | float, outer_radius_deg: int | float, ra0_deg: int | float, dec0_deg: int | float, indices: Literal[True]) Tuple[SkyModel, ndarray[Any, dtype[int64]]]

Filters the sky according to an inner and outer circle from the phase center

Parameters:
  • inner_radius_deg – Inner radius in degrees

  • outer_radius_deg – Outer radius in degrees

  • ra0_deg – Phase center right ascension

  • dec0_deg – Phase center declination

  • indices – Optional parameter, if set to True,

we also return the indices of the filtered sky copy :return sky: Filtered copy of the sky

filter_by_radius_euclidean_flat_approximation(inner_radius_deg: int | float, outer_radius_deg: int | float, ra0_deg: int | float, dec0_deg: int | float, indices: Literal[False] = False) SkyModel
filter_by_radius_euclidean_flat_approximation(inner_radius_deg: int | float, outer_radius_deg: int | float, ra0_deg: int | float, dec0_deg: int | float, indices: Literal[True]) Tuple[SkyModel, ndarray[Any, dtype[int64]]]

Filters sources within an annular region using a flat Euclidean distance approximation suitable for large datasets managed by Xarray.

This function is designed for scenarios where the dataset size precludes in-memory spherical geometry calculations. By leveraging a flat Euclidean approximation, it permits the use of Xarray’s out-of-core computation capabilities, thus bypassing the limitations imposed by the incompatibility of astropy.visualization.wcsaxes.SphericalCircle with Xarray’s data structures. Although this method trades off geometric accuracy against computational efficiency, it remains a practical choice for large angular fields where the curvature of the celestial sphere can be reasonably neglected.

Parameters

inner_radius_degIntFloat

The inner radius of the annular search region, in degrees.

outer_radius_degIntFloat

The outer radius of the annular search region, in degrees.

ra0_degIntFloat

The right ascension of the search region’s center, in degrees.

dec0_degIntFloat

The declination of the search region’s center, in degrees.

indicesbool, optional

If True, returns the indices of the filtered sources in addition to the SkyModel object. Defaults to False.

Returns

SkyModel or tuple of (SkyModel, NDArray[np.int_])

The filtered SkyModel object, and optionally the indices of the filtered sources if indices is set to True.

Raises

KaraboSkyModelError

If the sources attribute is not populated in the SkyModel instance prior to invoking this function.

Notes

Use this function for large sky models where a full spherical geometry calculation is not feasible due to memory constraints. It is particularly beneficial when working with Xarray and Dask, facilitating scalable data analysis on datasets that are too large to fit into memory.

classmethod get_GLEAM_Sky(min_freq: float | None = None, max_freq: float | None = None) _TSkyModel

Creates a SkyModel containing GLEAM sources from https://vizier.cfa.harvard.edu/ service with more than 300,000 sources on different frequency-bands.

The survey was created using MWA-telescope and covers the entire sky south of dec +30.

GLEAM’s frequency-range: 72-231 MHz

The required .fits file will get downloaded and cached on disk.

Args:

min_freq: Set min-frequency in Hz for pre-filtering. max_freq: Set max-frequency in Hz for pre-filtering.

Returns:

GLEAM sky as SkyModel.

classmethod get_MALS_DR1V3_Sky(min_freq: float | None = None, max_freq: float | None = None) _TSkyModel

Creates a SkyModel containing “MALS V3” catalogue, of 715,760 sources, where the catalogue and it’s information are from ‘https://mals.iucaa.in/’.

The MeerKAT Absorption Line Survey (MALS) consists of 1655 hrs of MeerKAT time (anticipated raw data ~ 1.7 PB) to carry out the most sensitive search of HI and OH absorption lines at 0<z<2, the redshift range over which most of the cosmic evolution in the star formation rate density takes place. The MALS survey is described in ‘Gupta et al. (2016)’.

MALS sky-coverage consists of 392 sources trackings between all RA-angles and DEC[-78.80,32.35] and different radius.

MALS’s frequency-range: 902-1644 MHz.

For publications, please honor their work by citing them as follows: - If you describe MALS or associated science, please cite ‘Gupta et al. 2016’. - If you use DR1 data products, please cite ‘Deka et al. 2024’.

Args:

min_freq: Set min-frequency in Hz for pre-filtering. max_freq: Set max-frequency in Hz for pre-filtering.

Returns:

MALS sky as SkyModel.

classmethod get_MIGHTEE_Sky(min_freq: float | None = None, max_freq: float | None = None) _TSkyModel

Creates a SkyModel containing “MIGHTEE Continuum Early Science L1 catalogue” sources from https://archive.sarao.ac.za service, consisting of 9896 sources.

MIGHTEE is an extragalactic radio-survey carried out using MeerKAT-telescope. It covers roughly the area between RA[149.40,150.84] and DEC[1.49,2.93].

MIGHTEES’s frequency-range: 1304-1375 MHz

The required .fits file will get downloaded and cached on disk.

Args:

min_freq: Set min-frequency in Hz for pre-filtering. max_freq: Set max-frequency in Hz for pre-filtering.

Returns:

MIGHTEE sky as SkyModel.

static get_OSKAR_sky(sky: ndarray[Any, dtype[float64]] | ndarray[Any, dtype[object_]] | DataArray | SkyModel, precision: Literal['single', 'double'] = 'double') Sky

Get OSKAR sky model object from the defined Sky Model

Returns:

oskar sky model

static get_fits_catalog(path: str) Table

Gets astropy.table.table.Table from a .fits catalog

Parameters:

path – Location of the .fits file

Returns:

fits catalog

classmethod get_sample_simulated_catalog() _TSkyModel

Downloads a sample simulated HI source catalog and generates a sky model using the downloaded data. The catalog size is around 8MB.

Source: The simulated catalog data was provided by Luis Machado (https://github.com/lmachadopolettivalle) in collaboration with the ETHZ Cosmology Research Group.

Returns:

SkyModel: The corresponding sky model. The sky model contains the following information:

  • ‘Right Ascension’ (ra): The right ascension coordinates

    of the celestial objects.

  • ‘Declination’ (dec): The declination coordinates of the

    celestial objects.

  • ‘Flux’ (i): The flux measurements of the celestial objects.

  • ‘Observed Redshift’: Additional observed redshift information

    of the celestial objects.

Note: Other properties such as ‘stokes_q’, ‘stokes_u’, ‘stokes_v’,

‘ref_freq’, ‘spectral_index’, ‘rm’, ‘major’, ‘minor’, ‘pa’, and ‘id’

are not included in the sky model.

classmethod get_sky_model_from_fits(*, fits_file: Path | str, prefix_mapping: SkyPrefixMapping, unit_mapping: Dict[str, UnitBase], units_sources: SkySourcesUnits | None = None, min_freq: float | None = None, max_freq: float | None = None, encoded_freq: UnitBase | None = None, chunks: int | Literal['auto'] | Tuple[int, ...] | Mapping[Hashable, int] | None = 'auto', memmap: bool = False) _TSkyModel

Creates a sky-model from fits_file. The following formats are supported: - Each data-array of the .fits file maps to a single SkyModel.sources column. - Frequency of some columns is encoded in col-names of the .fits file.

Args:

fits_file: The .fits file to create the sky-model from. prefix_mapping: Formattable col-names of .fits file. If encoded_freq is

not None, the freq-encoded field-values must have ‘{0}’ as placeholder.

unit_mapping: Mapping from col-unit to astropy.unit. units_sources: Units of SkyModel.sources. min_freq: Filter by min-freq in Hz? May increase file-reading significantly. max_freq: Filter by max-freq in Hz? May increase file-reading significantly. encoded_freq: Unit of col-name encoded frequency, if the .fits file has

it’s ref-frequency encoded in the col-names.

chunks: Coerce the array’s data into dask arrays with the given chunks.

Might be useful for reading larger-than-memory .fits files.

memmap: Whether to use memory mapping when opening the FITS file.

Allows for reading of larger-than-memory files.

Returns:

sky-model with according sources from fits_file.

classmethod get_sky_model_from_h5_to_xarray(path: str, prefix_mapping: SkyPrefixMapping | None = None, load_as: Literal['numpy_array', 'dask_array'] = 'dask_array', chunksize: int | Literal['auto'] = 'auto') _TSkyModel

Load a sky model dataset from an HDF5 file and converts it to an xarray DataArray.

Parameters

pathstr

Path to the input HDF5 file.

prefix_mappingSkyPrefixMapping

Mapping column names to their corresponding dataset paths in the HDF5 file. If the column is not present in the HDF5 file, set its value to None.

load_asLiteral[“numpy_array”, “dask_array”], default=”dask_array”

What type of array to load the data inside the xarray Data Array as.

chunksizeUnion[int, str], default=auto

Chunk size for Dask arrays. This determines the size of chunks that the data will be divided into when read from the file. Can be an integer or ‘auto’. If ‘auto’, Dask will choose an optimal chunk size.

Returns

xr.DataArray

A 2D xarray DataArray containing the sky model data. Rows represent data points and columns represent different data fields (‘ra’, ‘dec’, …).

get_wcs() WCS

Gets the currently active world coordinate system astropy.wcs For details see https://docs.astropy.org/en/stable/wcs/index.html

Returns:

world coordinate system

classmethod read_from_file(path: str) _TSkyModel

Read a CSV file in to create a SkyModel. The CSV should have the following columns

  • right ascension (deg)

  • declination (deg)

  • stokes I Flux (Jy)

  • stokes Q Flux (Jy): if no information available, set to 0

  • stokes U Flux (Jy): if no information available, set to 0

  • stokes V Flux (Jy): if no information available, set to 0

  • reference_frequency (Hz): if no information available, set to 0

  • spectral index (N/A): if no information available, set to 0

  • rotation measure (rad / m^2): if no information available, set to 0

  • major axis FWHM (arcsec): if no information available, set to 0

  • minor axis FWHM (arcsec): if no information available, set to 0

  • position angle (deg): if no information available, set to 0

  • true redshift: defaults to 0

  • observed redshift: defaults to 0

  • source id (object): is in SkyModel.source_ids if provided

Parameters:

path – file to read in

Returns:

SkyModel

static read_healpix_file_to_sky_model_array(file: str, channel: int, polarisation: Polarisation) Tuple[DataArray, int]

Read a healpix file in hdf5 format. The file should have the map keywords:

Parameters:
  • file – hdf5 file path (healpix format)

  • channel – Channels of observation (between 0 and maximum numbers of channels of observation)

  • polarisation – 0 = Stokes I, 1 = Stokes Q, 2 = Stokes U, 3 = Stokes V

Returns:

save_sky_model_as_csv(path: str) None

Save source array into a csv. :param path: path to save the csv file in.

set_wcs(wcs: WCS) None

Sets a new world coordinate system astropy.wcs For details see https://docs.astropy.org/en/stable/wcs/index.html

Parameters:

wcs – world coordinate system

setup_default_wcs(phase_center: List[int] | List[float] = [0, 0]) WCS

Defines a default world coordinate system astropy.wcs For more details see https://docs.astropy.org/en/stable/wcs/index.html

Parameters:

phase_center – ra-dec location

Returns:

wcs

classmethod sky_test() _TSkyModel

Construction of a sky model which can be used for testing and visualizing the simulation with equal distributed point sources around the phase center ra=20, deg=-30.

Returns:

The test sky model.

to_np_array(with_obj_ids: bool = False) ndarray[Any, dtype[float64]] | ndarray[Any, dtype[object_]]

Gets the sources as np.ndarray

Parameters:

with_obj_ids – Option whether object ids should be included or not

Returns:

the sources of the SkyModel as np.ndarray

to_sky_xarray(sources: ndarray[Any, dtype[float64]] | ndarray[Any, dtype[object_]] | DataArray) DataArray
Converts a np.ndarray or xr.DataArray to SkyModel.sources

compatible xr.DataArray.

Args:

sources: Array to convert. Col-order see SkyModel description.

Returns:

SkyModel.sources compatible xr.DataArray.

class BeamPattern(cst_file_path: str, telescope: Telescope | None = None, freq_hz: float = 0.0, pol: Literal['x', 'Y', 'XY'] = 'XY', element_type_index: int = 0, average_fractional_error_factor_increase: float = 1.1, ignore_data_at_pole: bool = True, avg_frac_error: float = 0.8, beam_method: Literal['Gaussian Beam', 'EIDOS_AH', 'EIDOS_EM', 'KatBeam'] = 'Gaussian Beam', interpol: Literal['inter2d', 'RectBivariateSpline'] = 'RectBivariateSpline')

:param

__init__(cst_file_path: str, telescope: Telescope | None = None, freq_hz: float = 0.0, pol: Literal['x', 'Y', 'XY'] = 'XY', element_type_index: int = 0, average_fractional_error_factor_increase: float = 1.1, ignore_data_at_pole: bool = True, avg_frac_error: float = 0.8, beam_method: Literal['Gaussian Beam', 'EIDOS_AH', 'EIDOS_EM', 'KatBeam'] = 'Gaussian Beam', interpol: Literal['inter2d', 'RectBivariateSpline'] = 'RectBivariateSpline') None
static get_eidos_holographic_beam(npix: int, ch: int, dia: int, thres: int, mode: Literal['AH', 'EM'] = 'AH') ndarray[Any, dtype[complex128]]

Returns beam

static get_meerkat_uhfbeam(f: int | float, pol: Literal['H', 'V', 'I'], beamextentx: int | float, beamextenty: int | float, sampling_step: int = 80) Tuple[ndarray[Any, dtype[float64]], ndarray[Any, dtype[float64]], ndarray[Any, dtype[float64]]]
Parameters:
  • pol

  • beamextent

Returns:

make_cst_from_arr(arr: _SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | bool | int | float | complex | str | bytes | _NestedSequence[bool | int | float | complex | str | bytes], output_file_path: str) None

Takes array of dimensions (*,8), and returns a cst files :param arr: :return: cst file with given output filename

plot_beam(theta: ndarray[Any, dtype[float64]], phi: ndarray[Any, dtype[float64]], absdir: ndarray[Any, dtype[float64]], path: str | None = None) None
Parameters:
  • theta – in radians

  • phi – in radian

  • absdir – in DBs

Returns:

polar plot

save_meerkat_cst_file(cstdata: ndarray[Any, dtype[float64]]) None

Save CST file for MeerKat telescope for the custom beams

static show_kat_beam(beampixels: ndarray[Any, dtype[float64]], beamextent: int, freq: int, pol: str, path: str | None = None) None
Parameters:
  • beamextent

  • freq

  • pol

Returns:

sim_beam(beam_method: Literal['Gaussian Beam', 'EIDOS_AH', 'EIDOS_EM', 'KatBeam'] | None = None, f: float | None = None, fov: int | float = 30.0, interpol: Literal['inter2d', 'RectBivariateSpline'] | None = None) Tuple[List[Quantity], ndarray[Any, dtype[float64]], ndarray[Any, dtype[float64]], ndarray[Any, dtype[float64]], ndarray[Any, dtype[float64]]]

Simulates the primary beam

Parameters:
  • beam_method – you can choose as beams: “Gaussian Beam”, “Eidos_AH”, “Eidos_EM”, “KatBeam”

  • f – the frequency for which the beam is simulated (MHz)

  • fov – when using “KatBeam” the beam will be simulated for this size of fov (radius)

class Observation(*, mode: str = 'Tracking', start_frequency_hz: int | float = 0, start_date_and_time: datetime | str, length: timedelta = datetime.timedelta(seconds=14400), number_of_channels: int = 1, frequency_increment_hz: int | float = 0, phase_centre_ra_deg: int | float = 0, phase_centre_dec_deg: int | float = 0, number_of_time_steps: int = 1)
class ObservationLong(*, mode: str = 'Tracking', start_frequency_hz: int | float = 0, start_date_and_time: datetime | str, length: timedelta = datetime.timedelta(seconds=14400), number_of_channels: int = 1, frequency_increment_hz: int | float = 0, phase_centre_ra_deg: int | float = 0, phase_centre_dec_deg: int | float = 0, number_of_time_steps: int = 1, number_of_days: int = 2)

This class allows the use of several observations on different days over a certain period of time within one day. If only ONE observation is desired, even if it takes a little longer, this is already possible using Observation. This class extends Observation so its parameters (except length) are not discussed here. length is little different, which describes the duration of ONE observation, whose maximum duration for ObservationLong is 24h.

Variables:

number_of_days – Number of successive days to observe

__init__(*, mode: str = 'Tracking', start_frequency_hz: int | float = 0, start_date_and_time: datetime | str, length: timedelta = datetime.timedelta(seconds=14400), number_of_channels: int = 1, frequency_increment_hz: int | float = 0, phase_centre_ra_deg: int | float = 0, phase_centre_dec_deg: int | float = 0, number_of_time_steps: int = 1, number_of_days: int = 2) None
Args:

start_date_and_time (Union[datetime, str]): Start time UTC and date for the observation. Strings are converted to datetime objects using datetime.fromisoformat.

mode (str, optional): TODO. Defaults to “Tracking”.

start_frequency_hz (IntFloat, optional): The frequency at the start of the first channel in Hz. Defaults to 0.

length (timedelta, optional): Length of observation. Defaults to timedelta(hours=4).

number_of_channels (int, optional): Number of channels / bands to use. Defaults to 1.

frequency_increment_hz (IntFloat, optional): Frequency increment between successive channels in Hz. Defaults to 0.

phase_centre_ra_deg (IntFloat, optional): Right Ascension of the observation pointing (phase centre) in degrees. Defaults to 0.

phase_centre_dec_deg (IntFloat, optional): Declination of the observation pointing (phase centre) in degrees. Defaults to 0.

number_of_time_steps (int, optional): Number of time steps in the output data during the observation length. This corresponds to the number of correlator dumps for interferometer simulations, and the number of beam pattern snapshots for beam pattern simulations. Defaults to 1.

class Telescope(longitude: float, latitude: float, altitude: float = 0.0)

WGS84 longitude and latitude and altitude in metres centre of the telescope.png centre. A telescope is described as follows:

Each row represents one station, with the elements being the horizontal x (east), horizontal y (north), and horizontal z (up) coordinates, followed by the errors in horizontal y (east), horizontal y (north), and horizontal z (up). Example: [[x, y, z, error_x, error_y, error_z], […]]

centre_longitudefloat

WGS84 longitude at the center of the telescope.

centre_latitudefloat

WGS84 latitude at the center of the telescope.

centre_altitudefloat

Altitude (in meters) at the center of the telescope.

__init__(longitude: float, latitude: float, altitude: float = 0.0) None

__init__ method

Parameters

longitudefloat

WGS84 longitude at the center of the telescope.

latitudefloat

WGS84 latitude at the center of the telescope.

altitudefloat, optional

Altitude (in meters) at the center of the telescope, default is 0.

add_antenna_to_station(station_index: int, horizontal_x: float, horizontal_y: float, horizontal_z: float = 0, horizontal_x_coordinate_error: float = 0, horizontal_y_coordinate_error: float = 0, horizontal_z_coordinate_error: float = 0) None

Add a new antenna to an existing station

Parameters:
  • station_index – Index of station to add antenna to

  • horizontal_x – east coordinate relative to the station center in metres

  • horizontal_y – north coordinate relative to the station center in metres

  • horizontal_z – altitude of antenna

  • horizontal_x_coordinate_error – east coordinate error

relative to the station center in metres :param horizontal_y_coordinate_error: north coordinate error relative to the station center in metres :param horizontal_z_coordinate_error: altitude of antenna error :return:

add_station(horizontal_x: float, horizontal_y: float, horizontal_z: float = 0.0, horizontal_x_coordinate_error: float = 0.0, horizontal_y_coordinate_error: float = 0.0, horizontal_z_coordinate_error: float = 0.0) None

Specify the stations as relative to the centre position :param horizontal_x: east coordinate relative to centre :param horizontal_y: north coordinate relative to centre :param horizontal_z: up coordinate :param horizontal_x_coordinate_error: east coordinate error :param horizontal_y_coordinate_error: north coordinate error :param horizontal_z_coordinate_error: up coordinate error

classmethod ang_res(freq: float, b: float) float

Calculates angular resolution in arcsec.

Angular resolution: θ=λ/B Wavelength: λ=c/f B: Max baseline in meters

Args:

freq: Frequency [Hz]. b: Max baseline in meters (e.g. from longest_baseline).

Returns:

Angular resolution in arcsec.

classmethod constructor(name: Literal['ACA', 'ALMA', 'ATCA', 'CARMA', 'NGVLA', 'PDBI', 'SMA', 'VLA'], version: Enum, backend: Literal[SimulatorBackend.OSKAR] = SimulatorBackend.OSKAR) Telescope
classmethod constructor(name: Literal['EXAMPLE', 'MeerKAT', 'ASKAP', 'LOFAR', 'MKATPlus', 'SKA1LOW', 'SKA1MID', 'VLBA', 'WSRT'], version: Literal[None] = None, backend: Literal[SimulatorBackend.OSKAR] = SimulatorBackend.OSKAR) Telescope
classmethod constructor(name: Literal['LOWBD2', 'LOWBD2-CORE', 'LOW', 'LOWR3', 'LOWR4', 'LOW-AA0.5', 'MID', 'MIDR5', 'MID-AA0.5', 'MEERKAT+', 'ASKAP', 'LOFAR', 'VLAA', 'VLAA_north'], version: Literal[None] = None, backend: Literal[SimulatorBackend.RASCIL] = SimulatorBackend.RASCIL) Telescope

Main constructor to obtain a pre-configured telescope instance. :param name: Name of the desired telescope configuration.

This name, together with the backend, is used as the key to look up the correct telescope specification file.

Parameters:
  • version – Version details required for some telescope configurations. Defaults to None.

  • backend – Underlying package to be used for the telescope configuration, since each package stores the arrays in a different format. Defaults to OSKAR.

Raises:

ValueError if the combination of input parameters is invalid. Specifically, if the requested telescope requires a version, but an invalid version (or no version) is provided, or if the requested telescope name is not supported by the requested backend.

Returns:

Telescope instance.

classmethod create_baseline_cut_telescope(lcut: bool | int | integer | float | floating, hcut: bool | int | integer | float | floating, tel: Telescope, tm_path: Path | str | None = None) Tuple[Path | str, Dict[str, str]]

Cut telescope tel for baseline-lengths.

Args:

lcut: Lower cut hcut: Higher cut tel: Telescope to cut off tm_path: .tm file-path to save the cut-telescope.

tm_path will get overwritten if it already exists.

Returns:

.tm file-path & station-name conversion (e.g. station055 -> station009)

get_OSKAR_telescope() Telescope

Retrieve the OSKAR Telescope object from the karabo.Telescope object.

Returns:

OSKAR Telescope object

get_baselines_dists() ndarray[Any, dtype[float64]]

Gets the interferometer baselines distances in meters.

It’s euclidean distance, not geodesic.

Returns:

Interferometer baselines dists in meters.

get_baselines_wgs84() ndarray[Any, dtype[float64]]

Gets the interferometer baselines in WGS84.

This function assumes that self.stations provides WGS84 coordinates.

Returns:

Baselines lon[deg]/lat[deg]/alt[m] (nx3).

longest_baseline() float64

Gets the longest baseline in meters.

Returns:

Length of longest baseline.

property name: str | None

Gets the telescope name (if available).

It’s just the file-name of the referred telescope-file without the ending.

Returns:

Telescope name or None.

plot_telescope(file: str | None = None) None

Plot the telescope according to which backend is being used, and save the resulting image into a file, if any is provided.

plot_telescope_OSKAR(file: str | None = None, block: bool = False) None

Plot the telescope and all its stations and antennas with longitude altitude

write_to_disk(dir_name: Path | str, *, overwrite: bool = False) None

Write dir_path to disk (must have .tm ending).

Parameters:
  • dir – directory in which the configuration will be saved in.

  • overwrite – If True an existing directory is overwritten if exists. Be careful to put the correct dir as input because the old one can get removed!

class InterferometerSimulation(ms_file_path: str | None = None, vis_path: str | None = None, channel_bandwidth_hz: int | float = 0, time_average_sec: int | float = 0, max_time_per_samples: int = 8, correlation_type: CorrelationType = CorrelationType.Cross_Correlations, uv_filter_min: int | float = 0, uv_filter_max: int | float = inf, uv_filter_units: FilterUnits = FilterUnits.WaveLengths, force_polarised_ms: bool = False, ignore_w_components: bool = False, noise_enable: bool = False, noise_seed: str | int = 'time', noise_start_freq: int | float = 1000000000.0, noise_inc_freq: int | float = 100000000.0, noise_number_freq: int = 24, noise_rms_start: int | float = 0, noise_rms_end: int | float = 0, noise_rms: Literal['Range', 'Data file', 'Telescope model'] = 'Range', noise_freq: Literal['Range', 'Data file', 'Observation settings', 'Telescope model'] = 'Range', enable_array_beam: bool = False, enable_numerical_beam: bool = False, beam_polX: BeamPattern | None = None, beam_polY: BeamPattern | None = None, use_gpus: bool | None = None, use_dask: bool | None = None, split_observation_by_channels: bool = False, n_split_channels: int | str = 'each', client: Client | None = None, precision: Literal['single', 'double'] = 'single', station_type: Literal['Aperture array', 'Isotropic beam', 'Gaussian beam', 'VLA (PBCOR)'] = 'Isotropic beam', enable_power_pattern: bool = False, gauss_beam_fwhm_deg: int | float = 0.0, gauss_ref_freq_hz: int | float = 0.0, ionosphere_fits_path: str | None = None, ionosphere_screen_type: str | None = None, ionosphere_screen_height_km: float | None = 300, ionosphere_screen_pixel_size_m: float | None = 0, ionosphere_isoplanatic_screen: bool | None = False)

Class containing all configuration for the Interferometer Simulation. :ivar ms_file_path: Optional. File path of the MS (mass spectrometry) file. :ivar vis_path: Optional. File path for visibilities output. :ivar channel_bandwidth_hz: The channel width, in Hz, used to simulate bandwidth

smearing. (Note that this can be different to the frequency increment if channels do not cover a contiguous frequency range.)

Variables:
  • time_average_sec – The correlator time-average duration, in seconds, used to simulate time averaging smearing.

  • max_time_per_samples – The maximum number of time samples held in memory before being written to disk.

  • correlation_type – The type of correlations to produce. Any value of Enum CorrelationType

  • uv_filter_min – The minimum value of the baseline UV length allowed by the filter. Values outside this range are not evaluated

  • uv_filter_max – The maximum value of the baseline UV length allowed by the filter. Values outside this range are not evaluated.

  • uv_filter_units – The units of the baseline UV length filter values. Any value of Enum FilterUnits

  • force_polarised_ms – If True, always write the Measurement Set in polarised format even if the simulation was run in the single polarisation ‘Scalar’ (or Stokes-I) mode. If False, the size of the polarisation dimension in the Measurement Set will be determined by the simulation mode.

  • ignore_w_components – If enabled, baseline W-coordinate component values will be set to 0. This will disable W-smearing. Use only if you know what you’re doing!

  • noise_enable – If true, noise is added.

  • noise_seed – Random number generator seed.

  • noise_start_freq – The start frequency in Hz for which noise is included, if noise is set to true.

  • noise_inc_freq – The frequency increment in Hz, if noise is set to true.

  • noise_number_freq – The number of frequency taken into account, if noise is set to true.

  • noise_rms_start – Station RMS (noise) flux density range start value, in Jy. The range is expanded linearly over the number of frequencies for which noise is defined.

  • noise_rms_end – Station RMS (noise) flux density range end value, in Jy. The range is expanded linearly over the number of frequencies for which noise is defined.

  • noise_rms

    The specifications for the RMS noise value:
    Telescope model: values are loaded from files in the telescope

    model directory.

    Data file: values are loaded from the specified file. Range: values are evaluated according to the specified range

    parameters (Default).

    The noise values are specified in Jy and represent the RMS noise of an unpolarised source in terms of flux measured in a single polarisation of the detector.

  • noise_freq

    The list of frequencies for which noise values are defined: Telescope model: frequencies are loaded from a data file in

    the telescope model directory.

    Observation settings: frequencies are defined by the observation

    settings.

    Data file: frequencies are loaded from the specified data file. Range: frequencies are specified by the range parameters

    (Default).

  • enable_array_beam – If true, then the contribution to the station beam from the array pattern (given by beam-forming the antennas in the station) is evaluated.

  • enable_numerical_beam – If true, make use of any available numerical element pattern files. If numerical pattern data are missing, the functional type will be used instead.

  • beam_polX – currently only considered for ObservationLong

  • beam_polX – currently only considered for ObservationLong

  • use_gpus – Set to true if you want to use gpus for the simulation

  • client – The dask client to use for the simulation

  • split_idxs_per_group – The indices of the sky model to split for each group of workers. If None, the sky model will not be split. Useful if the sky model is too large to fit into the memory of a single worker. Group index should be strictly monotonic increasing.

  • precision – For the arithmetic use you can choose between “single” or “double” precision

  • station_type – Here you can choose the type of each station in the interferometer. You can either disable all station beam effects by choosing “Isotropic beam”. Or select one of the following beam types: “Gaussian beam”, “Aperture array” or “VLA (PBCOR)”

  • enable_power_pattern – If true, gauss_beam_fwhm_deg will be taken in as power pattern.

  • gauss_beam_fwhm_deg – If you choose “Gaussian beam” as station type you need specify the full-width half maximum value at the reference frequency of the Gaussian beam here. Units = degrees. If enable_power_pattern is True, gauss_beam_fwhm_deg is in power pattern, otherwise it is in field pattern.

  • gauss_ref_freq_hz – The reference frequency of the specified FWHM, in Hz.

  • ionosphere_fits_path – The path to a fits file containing an ionospheric screen generated with ARatmospy. The file parameters (times/frequencies) should coincide with the planned observation.

__init__(ms_file_path: str | None = None, vis_path: str | None = None, channel_bandwidth_hz: int | float = 0, time_average_sec: int | float = 0, max_time_per_samples: int = 8, correlation_type: CorrelationType = CorrelationType.Cross_Correlations, uv_filter_min: int | float = 0, uv_filter_max: int | float = inf, uv_filter_units: FilterUnits = FilterUnits.WaveLengths, force_polarised_ms: bool = False, ignore_w_components: bool = False, noise_enable: bool = False, noise_seed: str | int = 'time', noise_start_freq: int | float = 1000000000.0, noise_inc_freq: int | float = 100000000.0, noise_number_freq: int = 24, noise_rms_start: int | float = 0, noise_rms_end: int | float = 0, noise_rms: Literal['Range', 'Data file', 'Telescope model'] = 'Range', noise_freq: Literal['Range', 'Data file', 'Observation settings', 'Telescope model'] = 'Range', enable_array_beam: bool = False, enable_numerical_beam: bool = False, beam_polX: BeamPattern | None = None, beam_polY: BeamPattern | None = None, use_gpus: bool | None = None, use_dask: bool | None = None, split_observation_by_channels: bool = False, n_split_channels: int | str = 'each', client: Client | None = None, precision: Literal['single', 'double'] = 'single', station_type: Literal['Aperture array', 'Isotropic beam', 'Gaussian beam', 'VLA (PBCOR)'] = 'Isotropic beam', enable_power_pattern: bool = False, gauss_beam_fwhm_deg: int | float = 0.0, gauss_ref_freq_hz: int | float = 0.0, ionosphere_fits_path: str | None = None, ionosphere_screen_type: str | None = None, ionosphere_screen_height_km: float | None = 300, ionosphere_screen_pixel_size_m: float | None = 0, ionosphere_isoplanatic_screen: bool | None = False) None
run_simulation(telescope: Telescope, sky: SkyModel, observation: Observation, backend: Literal[SimulatorBackend.OSKAR] = SimulatorBackend.OSKAR, primary_beam: None = None) Visibility
run_simulation(telescope: Telescope, sky: SkyModel, observation: ObservationLong, backend: Literal[SimulatorBackend.OSKAR] = SimulatorBackend.OSKAR, primary_beam: None = None) Visibility
run_simulation(telescope: Telescope, sky: SkyModel, observation: ObservationParallelized, backend: Literal[SimulatorBackend.OSKAR] = SimulatorBackend.OSKAR, primary_beam: None = None) List[Visibility]
run_simulation(telescope: Telescope, sky: SkyModel, observation: Observation, backend: Literal[SimulatorBackend.RASCIL], primary_beam: Image | None) Visibility
run_simulation(telescope: Telescope, sky: SkyModel, observation: ObservationAbstract, backend: SimulatorBackend, primary_beam: None = None) Visibility | List[Visibility] | Visibility

Run a single interferometer simulation with the given sky, telescope and observation settings. :param telescope: telescope model defining the configuration :param sky: sky model defining the sky sources :param observation: observation settings :param backend: Backend used to perform calculations (e.g. OSKAR, RASCIL) :param primary_beam: Primary beam to be included into visibilities.

Currently only relevant for RASCIL. For OSKAR, use the InterferometerSimulation constructor parameters instead.

set_ionosphere(file_path: str) None

Set the path to an ionosphere screen file generated with ARatmospy. The file parameters (times/frequencies) should coincide with the planned observation. see https://github.com/timcornwell/ARatmospy

Parameters:

file_path – file path to fits file.

simulate_foreground_vis(telescope: Telescope, foreground: SkyModel, foreground_observation: Observation, foreground_vis_file: str, write_ms: bool, foreground_ms_file: str) Tuple[Visibility, List[ndarray[Any, dtype[complex128]]], VisHeader, Binary, VisBlock, ndarray[Any, dtype[float64]], ndarray[Any, dtype[float64]], ndarray[Any, dtype[float64]]]

Simulates foreground sources