fastimgproto.telescope

Modules

fastimgproto.telescope.base module

class fastimgproto.telescope.base.Telescope(centre, ant_labels=None, ant_itrf_xyz=None, ant_local_xyz=None)[source]

Data-structure for holding together array layout and central position.

Also takes care of baseline calculation and antenna/baseline label-tracking.

Instantiate using the named constructors, e.g.:

meerkat = Telescope.from_itrf(
    *parse_itrf_ascii_to_xyz_and_labels(meerkat_itrf_filepath)
)
centre

astropy.coordinates.EarthLocation – Centre of the array

ant_labels

list – Antennae labels

ant_itrf_xyz

numpy.ndarray – Antennae xyz_positions in ITRF frame. [dtype: np.float_, shape (n_antennae,3,).

ant_local_xyz

numpy.ndarray – Antennae xyz_positions in local-XYZ frame. [dtype: np.float_, shape (n_antennae,3,).

baseline_local_xyz

numpy.ndarray – Baseline xyz-vectors in local-XYZ frame. [dtype: np.float_, shape (n_baselines,3,).

baseline_labels

list – Baseline labels. (See also generate_baselines_and_labels())

static from_itrf(ant_itrf_xyz, ant_labels)[source]

Instantiate telescope model from ITRF co-ordinates of antennae

(aka ECEF-coords, i.e. Earth-Centered-Earth-Fixed reference frame.)

Takes care of calculating central Latitude and Longitude from the mean antenna position, and converts the antenna positions into local-XYZ frame.

Parameters:
  • ant_itrf_xyz (numpy.ndarray) – Array co-ordinatates in the ITRF frame
  • ant_labels (list[str]) – Antennae labels
Returns:

A telescope class with the given array co-ords.

Return type:

Telescope

lha(ra, time)[source]

Calculate the local hour angle of a target-RA at a given time

Parameters:
Returns:

Local Hour Angle of target-RA.

Return type:

astropy.coordinates.Longitude

lst(time)[source]

Calculate the local sidereal time at the telescope

Parameters:time (astropy.time.Time) – Global timestamp
Returns:Local sidereal time expressed as an angle-Quantity.
Return type:astropy.coordinates.Longitude
next_transit(target_ra, start_time)[source]

Wrapper around time_of_next_transit()

See time_of_next_transit() for details.

Parameters:
Returns:

Approximate time of the next transit

Return type:

astropy.time.Time

uvw_at_local_hour_angle(local_hour_angle, dec)[source]

Transform baselines to UVW for source at given local-hour-angle and Dec.

Parameters:
  • local_hour_angle (astropy.units.Quantity) – Local hour angle. LHA=0 is transit, LHA=-6h is rising, LHA=+6h is setting. [Dimensions of angular width / arc, e.g. u.deg, u.rad, or u.hourangle]
  • dec (astropy.units.Quantity) – Declination. [Dimensions of angular width / arc]
Returns:

UVW-array, with units of metres. [Shape: (3, n_baselines), dtype: numpy.float_].

Return type:

astropy.units.Quantity

uvw_tracking_skycoord(pointing_centre, obs_times, progress_bar=None)[source]

Calculate the UVW-array towards pointing centre for all obs_times.

This function calculates an ndarray of uvw-coords corresponding to an instantaneous observation at each of obs_times, concatenating the results to produce an ndarray of shape (n_baselines * len(obs_times), 3)

Parameters:
  • pointing_centre (astropy.coordinates.SkyCoord) – Pointing centre co-ords for UVW calculation.
  • obs_times (list) – List of astropy.time.Time, the instants of observation.
  • progress_bar (tqdm.tqdm) – [Optional] progressbar to update.
Returns:

UVW-array, with units of metres.

Return type:

astropy.units.Quantity

fastimgproto.telescope.base.generate_baselines_and_labels(antenna_positions, antenna_labels)[source]

Given an array of antenna positions, generate an array of baseline-vectors.

Baseline-vectors are generated as follows:

  • It is assumed that the antenna-positions supplied are already in the desired order, with matching labels in the antenna_labels list. (Typically something like ['ANT1','ANT2','ANT3'].

  • Pairings are then generated according to position in the supplied list, as follows:

    pairings = list(itertools.combinations(range(len(antenna_positions)), 2))
    
  • Baseline labels are concatenated with a comma, e.g. pairing of 'ANT1' and 'ANT2' has label 'ANT1,ANT2'.

  • Baseline vectors represent translation from first antenna to second antenna, therefore:

    baseline_1_2 = ant_2_pos - ant_1_pos
    
Parameters:
  • antenna_positions (astropy.units.Quantity) – Antennae xyz_positions in local-XYZ frame. [dtype: np.float_, shape (n_antennae,3,).
  • antenna_labels (list) – Antennae labels
Returns:

(baseline_vecs, baseline_labels). Both of length \(n_{baselines} = {n_{antennae} \choose 2}\), where baseline_labels is simply a list of strings, and baseline_vecs is a (numpy.ndarray) [dtype: np.float_, shape (n_baselines,3,).

Return type:

tuple

fastimgproto.telescope.base.hstack_table_columns_as_ndarray(cols)[source]

When passed a subset of astropy Table columns, turn them into a 2-d ndarray.

(Note, columns should have same dtype)

Example

>>> import numpy as np
>>> import astropy.table
>>> x = np.arange(5)
>>> y = x**2
>>> t = astropy.table.Table(data=(x,y), names=('x','y'))
>>> hstack_table_columns_as_ndarray(t.columns[:2])
array([[ 0,  0],
       [ 1,  1],
       [ 2,  4],
       [ 3,  9],
       [ 4, 16]])
Parameters:cols (list) – List (or iterable) of astropy.table.TableColumns. Typically a table or a slice of a table’s columns.
Returns:H-stacked columns-array.
Return type:numpy.ndarray
fastimgproto.telescope.base.parse_itrf_ascii_to_xyz_and_labels(tabledata)[source]

Read ITRF data in ascii format and return (antenna positions, labels).

ASCII table should be simple whitespace-delimited, column order: x y z diameter label mount

Parameters:tabledata – (str, file-like, or list). File-object or path-to-file.
Returns:Tuple of (ant_itrf_xyz, ant_labels), where ant_itrf_xyz is a numpy.ndarray of shape (n,3), dtype np.float_ and ant_labels is a numpy.ndarray of shape (n,), dtype np.str_.
Return type:tuple