Array Telescope Simulation (telsim)¶
- simms.telescope.generate_ms.create_ms(ms: str, telescope_name: str | File, pointing_direction: str, dtime: int, ntimes: int, start_freq: str | float, dfreq: str | float, nchan: int, correlations: str, row_chunks: int, sefd: float, column: str, smooth: str = None, fit_order: int = None, start_time: str | List[str] = None, start_ha: float = None, freq_range: str = None, sfile: File = None, tsys_over_eta: float = None, subarray_list: List[str] = None, subarray_range: List[int] = None, subarray_file: File = None, low_source_limit: float | str = None, high_source_limit: float | str = None)¶
Generate a CASA Measurement Set (MS) for a simulated observation.
- Parameters:
ms (str) – Path to output Measurement Set directory (will be overwritten if it exists).
telescope_name (Union[str, File]) – Telescope or array definition (name or configuration file).
pointing_direction (str) – Target sky direction (e.g. ‘J2000 12:34:56.7 -45.23.10’).
dtime (int) – Integration time per visibility (seconds).
ntimes (int) – Number of time steps.
start_freq (Union[str, float]) – Starting frequency (string with units or raw Hz).
dfreq (Union[str, float]) – Channel frequency increment (string with units or raw Hz).
nchan (int) – Number of spectral channels.
correlations (str) – Correlation types (e.g. ‘XX’, ‘YY’, ‘RR’, ‘LL’).
row_chunks (int) – Dask chunk size along the row (time-baseline) dimension.
sefd (float) – System Equivalent Flux Density (Jy). If provided, noise is generated.
column (str) – Name of data column to populate (e.g. ‘DATA’ or ‘CORRECTED_DATA’).
smooth (str, optional) – SEFD frequency smoothing method (‘polyn’ or ‘spline’).
fit_order (int, optional) – Polynomial degree or spline smoothing parameter.
start_time (Union[str, List[str]], optional) – Observation start time (string or list for multi-part definitions).
start_ha (float, optional) – Starting hour angle (radians).
freq_range (str, optional) – Alternative frequency specification: [start, end, nchan].
sfile (File, optional) – Sensitivity file to override SEFD values.
tsys_over_eta (float, optional) – Tsys/eta value used when deriving SEFD if not explicitly given.
subarray_list (List[str], optional) – List of antenna names to include.
subarray_range (List[int], optional) – Index range of antennas to include.
subarray_file (File, optional) – File specifying subarray selection.
low_source_limit (Union[float, str], optional) – Elevation (deg) below which rows are flagged.
high_source_limit (Union[float, str], optional) – Elevation (deg) above which rows are flagged.
- Returns:
Writes the Measurement Set and its subtables to disk.
- Return type:
None
- simms.telescope.generate_ms.nda(items)¶
Convert a list of items to a Dask array with appropriate dtype.
- simms.telescope.generate_ms.parse_frequency(value, option: str)¶
Parses a frequency value which may be a string with units or a float without units.
- Parameters:
value (Union[str, float, int]) – The frequency value to parse.
- simms.telescope.generate_ms.remove_ms(ms: File | str)¶
- class simms.telescope.array_utilities.Array(layout: str | File, degrees: bool = True, sefd: int | float | List[int | float] = None, tsys_over_eta: int | float | List[int | float] = None, sensitivity_file: File = None, subarray_list: List[str] = None, subarray_range: List[int] = None, subarray_file: File = None)¶
Bases:
objectThe Array class has functions for converting from one coordinate system to another.
- geodetic2global()¶
Convert the antenna positions from the geodetic (WGS84) frame to the global (ITRF/ECEF) frame
- Return type:
An array of the antennas XYZ positions and an array of the array center in XYZ
- geodetic2local()¶
Converts the antenna positions from the geodetic (WGS84) frame to the local (ENU) frame
- Return type:
An array of the antenna positions in the local frame ENU
- get_antenna_elevation(ant_latitudes: List[float], dec: float, h0: float, ntimes: int) List[float]¶
Get the elevations of the antennas at the observing times.
- Parameters:
- Returns:
antenna_elevations – : the elevations of the antennas in degrees.
- Return type:
Array[float]
- get_source_elevation(latitude: float, declination: float, hour_angles: List[float]) List[float]¶
Track the source during the observation and get its elevation.
- Parameters:
latitude (float) – : latitude of the observer in radians.
declination (float) – : declination of the source in radians.
hour_angles (array[float]) – : hour angles of the source in radians.
- Returns:
source_elevations – : the elevations of the source in degrees.
- Return type:
array[float]
- global2geodetic(X, Y, Z)¶
Convert the antenna positions from the global (ITRF/ECEF) frame to the geodetic (WGS84) frame.
- Parameters:
(X (float) – : global coordinates of the antennas.
Y (float) – : global coordinates of the antennas.
Z) (float) – : global coordinates of the antennas.
- Returns:
(phi,lam,h) – : geodetic coordinates of the antennas.
- Return type:
float
- set_itrf()¶
- uvgen(pointing_direction, dtime, ntimes, start_freq, dfreq, nchan, start_time=None, start_ha=None) ObjDict¶
Generate uvw coordimates
- Parameters:
pointing_direction (List[str]) – : pointing direction.
dtime (int) – : integration time.
ntimes (int) – : number of times the sky should be snapped.
start_freq (Union[str,float]) – : starting freq of the observation.
dfreq (Union[str,float]) – : frequency interval.
nchan (int) – : number of channels.
start_time (Union[str, List[str]]) – : start time of the observation in UTC (YYYY-MM-DDTHH:MM:SS). Default is the current machine time.
start_ha (float) – : start hour angle in radians
- Return type:
An array of the uvw time dependent positions, the time array and the frequency array
- simms.telescope.array_utilities.calculate_array_centre(ant_locations: List[List[float]]) List[float]¶
Calculate the geometric centre of an array given the antenna locations.
- Parameters:
ant_locations (List[List[float]]) – A list of antenna locations, where each location is a list of [X, Y, Z] coordinates.
- Returns:
The [X, Y, Z] coordinates of the geometric centre of the array.
- Return type:
List[float]
- simms.telescope.array_utilities.ms_addrow(ms, subtable: str, nrows: int)¶
- simms.telescope.array_utilities.write_centre_to_array_config(array_config_path: File, centre: List[float]) None¶
Write the calculated centre to the array configuration YAML file.
- Parameters:
array_config_path (File) – Path to the array configuration YAML file.
centre (List[float]) – The [X, Y, Z] coordinates of the geometric centre to write.