coordio Reference

High Level Utility Functions

class coordio.utils.Moffat2dInterp(Noffset=None, FWHM=None, beta=None)[source]

Bases: object

Create the dict of 1D interpolations function for a moffat profile offset for various FWHMs. Object returned includes two dicts, one for APO and one for LCO for the various fiber sizes.

class coordio.utils.MoffatLossProfile(offset, beta, FWHM, rfiber=1)[source]

Bases: object

class for calculating fiber loss based on Moffat profile

Parameters:
  • offset (numpy.array) – fiber offset

  • beta (float) – Power index of the Moffat profile

  • FWHM (float) – seeing

  • rfiber (float) – radius of the fiber

flux_loss(offset)[source]

Function computing the flux loss obtained by moving the fiber across the source in a certain direction.

Prameters

offset: numpy.array

fiber offset

returns:

norm (numpy.array) – the flux loss

func_magloss()[source]

Function computing the magnitude loss obtained by moving the fiber across the source in a certain direction. Note that the magnitude loss with offset = 0 is not 0, and should coincide with the difference between PSF and fiber magnitudes.

Prameters

returns:

magloss (numpy.array) – the magnitude loss

moffat_1D(r, amplitude=None)[source]

Function computing Moffat 1D profile

Parameters:

r (float or numpy.array) – Distance from the centre of the profile

Returns:

moff_prof (float or numpy.array) – 1D Moffat profile.

moffat_norm(amplitude)[source]

Function computing the normalized the Moffat 1D profile.

Parameters:

amplitude (float) – Amplitude of the profile

Returns:

moff_prof_norm (float) – Normalised Moffat profile.

coordio.utils.fitsTableToPandas(recarray)[source]
coordio.utils.gaia_mags2sdss_gri(gaia_g, gaia_bp=None, gaia_rp=None, gaia_bp_rp=None)[source]

Stolen from Tom Dwelly’s dither_work:

https://github.com/sdss/dither_work/blob/main/python/dither_work/utils.py#L75

There are some fairly accurate+reliable G,BP-RP -> SDSS gri transforms available in the literature. Evans et al (2018, https://www.aanda.org/articles/aa/full_html/2018/08/aa32756-18/aa32756-18.html) give transforms for main sequence stars. Paraphrased below:

sdss_g_from_gdr2 = phot_g_mean_mag_gdr2 - (0.13518 - 0.46245*bp_rp_gdr2 +

-0.25171*bp_rp_gdr2**2 + 0.021349*bp_rp_gdr2**3)

sdss_r_from_gdr2 = phot_g_mean_mag_gdr2 - (-0.12879 + 0.24662*bp_rp_gdr2 +

-0.027464*bp_rp_gdr2**2 - 0.049465*bp_rp_gdr2**3)

sdss_i_from_gdr2 = phot_g_mean_mag_gdr2 - (-0.29676 + 0.64728*bp_rp_gdr2 +

-0.10141*bp_rp_gdr2**2)

coordio.utils.object_offset(mags, mag_limits, lunation, waveName, obsSite, fmagloss=None, safety_factor=None, beta=None, FWHM=None, skybrightness=None, offset_min_skybrightness=None, can_offset=None)[source]

Returns the offset needed for object with mag to be observed at mag_limit. Currently assumption is all offsets are set in positive RA direction.

Parameters:
  • mags (numpy.array) – The magniutdes of the objects. Should be a 2D, Nx10 array, where N is number of objects and length of 10 index should correspond to magntidues: [g, r, i, z, bp, gaia_g, rp, J, H, K].

  • mag_limits (numpy.array) – Magnitude limits for the designmode of the design. This should be an array of length N=10 where indexes correspond to magntidues: [g, r, i, z, bp, gaia_g, rp, J, H, K]. This matches the apogee_bright_limit_targets_min or boss_bright_limit_targets_min (depending on instrument) from targetdb.DesignMode for the design_mode of the design.

  • lunation (str:) – If the designmode is bright time (‘bright’) or dark time (‘dark’)

  • waveName (str) – Instrument for the fibers offset definition being applied to. Either ‘Boss’ or ‘Apogee’.

  • obsSite (str) – The observatory of the observation. Should either be ‘APO’ or ‘LCO’.

  • fmagloss (object) – Moffat2dInterp class with the lookup table for doing the moffat profile inversion. If None, then table is calculated at function call.

  • safety_factor (float) – Factor to add to mag_limit. If None, default set in code.

  • beta (float) – Power index of the Moffat profile. If None, default set in code.

  • FWHM (float) – seeing for the Moffat profile. If None, default set in code.

  • skybrightness (float) – Sky brightness for the field cadence. Only set if want to check for offset_flag TOO_DARK (8).

  • offset_min_skybrightness (float) – Minimum sky brightness for the offset. Only set if want to check for offset_flag TOO_DARK (8).

  • can_offset (boolean or numpy.array) – can_offset value from targetdb for the target(s) to be offset. Only set if want to check for offset_flag NO_CAN_OFFSET (16).

Returns:

  • delta_ra (numpy.array) – offset in RA in arcseconds around object(s) in numpy array of length N

  • delta_dec (numpy.array) – offset in Decl. in arcseconds around object(s) in numpy array of length N

  • offset_flag (numpy.array) – bitmask for how offset was set in numpy array of length N. Flags are:

    • 0: offset applied normally (i.e. when mag <= mag_limit)

    • 1: no offset applied because mag > mag_limit

    • 2: no offset applied because magnitude was null value.

    • 8: offsets should not be used as sky brightness is <=

      minimum offset sky brightness

    • 16: no offsets applied because can_offset = False

    • 32: no offset applied because mag <= offset_bright_limit
      (offset_bright_limit is G = 6 for Boss bright time and

      G = 13 for Boss dark time, and H = 1 for Apogee).

    • 64: no offset applied because no valid magnitude limits

coordio.utils.offset_definition(mag, mag_limits, lunation, waveName, obsSite, fmagloss=None, safety_factor=0.0, beta=None, FWHM=None, skybrightness=None, offset_min_skybrightness=None, can_offset=None, use_type='bright_neigh', mag_limit_ind=None)[source]

Returns the offset needed for object with mag to be observed at mag_limit. This is for the piecewise appromixation used based on: https://wiki.sdss.org/pages/viewpage.action?pageId=100173069

Parameters:
  • mag (float or numpy.array) – The magniutde(s) of the objects. For BOSS should be Gaia G-band and for APOGEE should be 2MASS H-band.

  • mag_limits (numpy.array) – Magnitude limits for the designmode of the design. This should be an array of length N=10 where indexes correspond to magntidues: [g, r, i, z, bp, gaia_g, rp, J, H, K]. This matches the apogee_bright_limit_targets_min or boss_bright_limit_targets_min (depending on instrument) from targetdb.DesignMode for the design_mode of the design.

  • lunation (str:) – If the designmode is bright time (‘bright’) or dark time (‘dark’)

  • waveName (str) – Instrument for the fibers offset definition being applied to. Either ‘Boss’ or ‘Apogee’.

  • obsSite (str) – The observatory of the observation. Should either be ‘APO’ or ‘LCO’.

  • fmagloss (object) – Moffat2dInterp class with the lookup table for doing the moffat profile inversion. If None, then table is calculated at function call.

  • safety_factor (float) – Factor to add to mag_limit. Should equal zero for bright neighbor checks (i.e. remain at default).

  • beta (float) – Power index of the Moffat profile. If None, default set in code.

  • FWHM (float) – seeing for the Moffat profile. If None, default set in code.

  • skybrightness (float) – Sky brightness for the field cadence. Only set if want to check for offset_flag TOO_DARK (8).

  • offset_min_skybrightness (float) – Minimum sky brightness for the offset. Only set if want to check for offset_flag TOO_DARK (8).

  • can_offset (boolean or numpy.array) – can_offset value from targetdb for the target(s) to be offset. Only set if want to check for offset_flag NO_CAN_OFFSET (16).

  • use_type (str) – Defines type for definition use. ‘bright_neigh’ is for bright neighbor check and auto sets mag_limit from mag_limits array. ‘offfset’ is for offsetting and the index from mag_limits array is set by mag_limit_ind.

  • mag_limit_ind (int) – when used with use_type=’offfset’, then this sets what index from mag_limits array is set as the magnitude limit. Otherwise, for bright neighbor check this is set in code.

Returns:

  • r (float or numpy.array) – offset radius in arcseconds around object(s)

  • offset_flag (int or numpy.array) –

    bitmask for how offset was set. Flags are:
    • 0: offset applied normally (i.e. when mag <= mag_limit)

    • 1: no offset applied because mag > mag_limit

    • 2: no offset applied because magnitude was null value.

    • 8: No offset applied because sky brightness is <=

      minimum offset sky brightness

    • 16: no offsets applied because can_offset = False

    • 32: no offset applied because mag <= offset_bright_limit
      (offset_bright_limit is G = 6 for Boss bright time and

      G = 13 for Boss dark time, and H = 1 for Apogee).

coordio.utils.radec2wokxy(ra, dec, coordEpoch, waveName, raCen, decCen, obsAngle, obsSite, obsTime, focalScale=None, pmra=None, pmdec=None, parallax=None, radVel=None, pressure=None, relativeHumidity=0.5, temperature=10, fullOutput=False, darLambda=None)[source]

Convert from ra, dec ICRS coords to a flat-wok XY in mm. At obsAngle=0 wok +y is a aligned with +dec, and wok +x is aligned with +ra

Question for José, do we need to enforce a time scale? I think everything defaults to UTC.

Parameters:
  • ra (numpy.ndarray) – Right ascension in ICRS, degrees

  • dec (numpy.ndarray) – Declination in ICRS, degrees

  • coordEpoch (float) – A TDB Julian date, the epoch for the input coordinates (from the catalog). Defaults to J2000.

  • waveName (str, or numpy.ndarray) – Array elements must be “Apogee”, “Boss”, or “GFA” strings

  • raCen (float) – Right ascension of field center, in degrees

  • decCen (float) – Declination of field center, in degrees.

  • obsAngle (float) – Position angle for observation. Angle is measured from North through East to wok +y. So obsAngle of 45 deg, wok +y points NE.

  • obsSite (str) – Either “APO” or “LCO”

  • obsTime (float) – TDB Julian date. The time at which these coordinates will be observed with the FPS.

  • focalScale (float or None) – Scale factor for conversion between focal and wok coords, found via dither analysis. Defaults to value in defaults.SITE_TO_SCALE

  • pmra (numpy.ndarray) – A 1D array with the proper motion in the RA axis for the N targets, in milliarcsec/yr. Must be a true angle, i.e, it must include the cos(dec) term. Defaults to 0.

  • pmdec (numpy.ndarray) – A 1D array with the proper motion in the RA axis for the N targets, in milliarcsec/yr. Defaults to 0.

  • parallax (numpy.ndarray) – A 1D array with the parallax for the N targets, in milliarcsec. Defaults to 0.

  • radVel (numpy.ndarray) – A 1D array with the radial velocity in km/s, positive when receding. Defaults to 0.

  • pressure (float) –

    The atmospheric pressure at the site, in millibar (same as hPa). If not provided the pressure will be calculated using the altitude \(h\) and the approximate expression

    \[p \sim -p_0 \exp\left( \dfrac{g h M} {T_0 R_0} \right)\]

    where \(p_0\) is the pressure at sea level, \(M\) is the molar mass of the air, \(R_0\) is the universal gas constant, and \(T_0=288.16\,{\rm K}\) is the standard sea-level temperature.

  • relativeHumidity (float) – The relative humidity, in the range \(0-1\). Defaults to 0.5.

  • temperature (float) – The site temperature, in degrees Celsius. Defaults to \(10^\circ{\rm C}\).

  • fullOutput (Bool) – If True all intermediate coords are returned

Returns:

  • xWok (numpy.ndarray) – x wok coordinate, mm

  • yWok (numpy.ndarray) – y wok coordinate, mm

  • fieldWarn (numpy.ndarray) – boolean array. Where True the coordinate converted should be eyed with suspicion. (It came from very far off axis).

  • hourAngle (float) – hour angle of field center in degrees

  • positionAngle (float) – position angle of field center in degrees

  • if fullOutput=True additionally returns

  • xFocal (numpy.ndarray) – x focal coordinate, mm

  • yFocal (numpy.ndarray) – y focal coordinate, mm

  • thetaField (numpy.ndarray) – azimuthal field coordinate, deg +RA through +Dec

  • phiField (numpy.ndarray) – polar field coordinate, deg off axis

  • altitude (numpy.ndarray) – altitude coordinate, deg

  • azimuth (numpy.ndarray) – azimuth coordinate N=0, E=90, deg

  • altFieldCen (float) – altitude of field center, deg

  • azFieldCen (float) – azimuth of field center, deg

coordio.utils.refinexy(image, x, y, psf_sigma=2.0, cutout=19)[source]

Refines positions of stars in an image. :param image: 2-D ndarray :type image: numpy.float32 :param x: 1-D ndarray of rough x positions :type x: numpy.float32 :param y: 1-D ndarray of rough y positions :type y: numpy.float32 :param psf_sigma: sigma of Gaussian PSF to assume (default 2 pixels) :type psf_sigma: float :param cutout: size of cutout used, should be odd (default 19) :type cutout: int

Returns:

  • xr (ndarray of numpy.float32) – refined x positions

  • yr (ndarray of numpy.float32) – refined y positions

Notes

Calls drefine.c in libdimage.so copied directly from: https://github.com/blanton144/dimage

coordio.utils.simplexy(image, psf_sigma=1.0, plim=8.0, dlim=1.0, saddle=3.0, maxper=1000, maxnpeaks=5000)[source]

Determines positions of stars in an image. :param image: 2-D ndarray :type image: numpy.float32 :param psf_sigma: sigma of Gaussian PSF to assume (default 1 pixel) :type psf_sigma: float :param plim: significance to select objects on (default 8) :type plim: float :param dlim: tolerance for closeness of pairs of objects (default 1 pixel) :type dlim: float :param saddle: tolerance for depth of saddle point to separate sources

(default 3 sigma)

Parameters:
  • maxper (int) – maximum number of children per parent (default 1000)

  • maxnpeaks (int) – maximum number of stars to find total (default 100000)

Returns:

(x, y, flux) ((numpy.float32, numpy.float32, numpy.float32)) – ndarrays with pixel positions and peak pixel values of stars

Notes

Calls simplexy.c in libdimage.so copied directly from: https://github.com/blanton144/dimage

coordio.utils.wokCurveAPO(r)[source]

Curve of the wok at APO at radial position r

Parameters:
  • r (scalar or 1D array) – radius (cylindrical coords) mm

  • Returns

  • ---------

  • result (scalar or 1D array) – z (height) of wok surface in mm (0 at vertex)

Curve that was specified for machining the APO wok, decided by Colby J.

coordio.utils.wokCurveLCO(r)[source]

Curve of the wok LCO at radial position r

Parameters:
  • r (scalar or 1D array) – radius (cylindrical coords) mm

  • Returns

  • ---------

  • result (scalar or 1D array) – z (height) of wok surface in mm (0 at vertex)

  • wok (Curve that was specified for machining the LCO)

  • J. (decided by Colby)

coordio.utils.wokxy2radec(xWok, yWok, waveName, raCen, decCen, obsAngle, obsSite, obsTime, focalScale=None, pressure=None, relativeHumidity=0.5, temperature=10)[source]

Convert from flat-wok XY (mm) to ra, dec ICRS coords (deg)

Question for José, do we need to enforce a time scale? I think everything defaults to UTC.

Parameters:
  • xWok (numpy.ndarray) – X wok coordinates, mm

  • yWok (numpy.ndarray) – Y wok coordinates, mm

  • waveName (str, or numpy.ndarray) – Array elements must be “Apogee”, “Boss”, or “GFA” strings

  • raCen (float) – Right ascension of field center, in degrees

  • decCen (float) – Declination of field center, in degrees.

  • obsAngle (float) – Position angle for observation. Angle is measured from North through East to wok +y. So obsAngle of 45 deg, wok +y points NE.

  • obsSite (str) – Either “APO” or “LCO”

  • obsTime (float) – TDB Julian date. The time at which these coordinates will be observed with the FPS.

  • focalScale (float or None) – Scale factor for conversion between focal and wok coords, found via dither analysis. Defaults to value in defaults.SITE_TO_SCALE

  • pressure (float) –

    The atmospheric pressure at the site, in millibar (same as hPa). If not provided the pressure will be calculated using the altitude \(h\) and the approximate expression

    \[p \sim -p_0 \exp\left( \dfrac{g h M} {T_0 R_0} \right)\]

    where \(p_0\) is the pressure at sea level, \(M\) is the molar mass of the air, \(R_0\) is the universal gas constant, and \(T_0=288.16\,{\rm K}\) is the standard sea-level temperature.

  • relativeHumidity (float) – The relative humidity, in the range \(0-1\). Defaults to 0.5.

  • temperature (float) – The site temperature, in degrees Celsius. Defaults to \(10^\circ{\rm C}\).

Returns:

  • xWok (numpy.ndarray) – x wok coordinate, mm

  • yWok (numpy.ndarray) – y wok coordinate, mm

  • fieldWarn (numpy.ndarray) – boolean array. Where True the coordinate converted should be eyed with suspicion. (It came from very far off axis).

SOFA

class coordio.sofa_bindings.SOFA(name=None)[source]

Bases: CDLL

Load the SOFA C extension.

Parameters:

name (str) – The path to the shared object. If None, will try the standard variants of libsofa for this architecture. If the library is not in a standard location the path needs to be added to LD_LIBRARY_PATHS.

get_internal_date(date=None, scale='UTC')[source]

Returns the internal representation of a date.

Parameters:
  • date (datetime) – The date to convert. If None, utcnow will be used and scale will be ignored.

  • scale (str) – The scale of the date.

Returns:

d1, d2 (float) – The two-part Julian date.

Site

class coordio.site.Site(name, latitude=None, longitude=None, altitude=None, pressure=None, temperature=None, rh=None, system_time_scale=None)[source]

Bases: object

A representation of an observing site.

Parameters:
  • name (str) – The name of the site. If the name matches that of a site in the configuration file those values will be used, overridden by the input keywords.

  • latitude (float) – The site latitude in degrees.

  • longitude (float) – The East-positive longitude of the site, in degrees.

  • altitude (float) – The altitude (elevation) of the site above see level, in meters. Defaults to zero meters.

  • pressure (float) –

    The atmospheric pressure at the site, in millibar (same as hPa). If not provided the pressure will be calculated using the altitude \(h\) and the approximate expression

    \[p \sim -p_0 \exp\left( \dfrac{g h M} {T_0 R_0} \right)\]

    where \(p_0\) is the pressure at sea level, \(M\) is the molar mass of the air, \(R_0\) is the universal gas constant, and \(T_0=288.16\,{\rm K}\) is the standard sea-level temperature.

  • temperature (float) – The site temperature, in degrees Celsius. Defaults to \(10^\circ{\rm C}\).

  • rh (float) – The relative humidity, in the range \(0-1\). Defaults to 0.5.

  • system_time_scale (str) – The time scale of the system time. Defaults to UTC.

Variables:

time (.Time) – A Time instance describing the time for the observation. Use set_time to set the time.

copy()[source]

Returns a copy of the object.

set_time(value=None, scale='TAI')[source]

Sets the time of the observation.

Parameters:
  • value (float, .Time, or None) – A Time instance or the Julian date of the time to use. If None, the current system datetime will be used.

  • scale (str) – The time scale of the Julian date. Ignored if value is a Time instance.

IERS

class coordio.iers.IERS(path=None, channel='finals', download=True)[source]

Bases: object

Wrapper around the IERS bulletins.

Parameters:
  • path (str) – The path where the IERS table, in CSV format, lives or will be saved to. If None, defaults to config['iers']['path'].

  • channel (str) – The IERS channel to use. Right now only finals is implemented.

  • download (bool) – If the IERS table is not found on disk or it’s out of date, download an updated copy.

Variables:

data (ndarray) – A record array with the parsed IERS data, trimmed to remove rows for which UT1-UTC is not defined.

get_delta_ut1_utc(jd=None, download=True)[source]

Returns the interpolated UT1-UTC value, in seconds.

is_valid(jd=None, offset=5, download=True)[source]

Determines whether a JD is included in the loaded IERS table.

Parameters:
  • jd (float) – The Julian date to check. There is a certain ambiguity in the format of the MJD in the IERS tables but for most purposes the difference between UTC and TAI should be meaningless for these purposes.

  • offset (int) – Number of days to look ahead. If jd+offset is not included in the IERS table the method will return False.

  • download (bool) – Whether to automatically download an updated version of the IERS table if the requested date is not included in the current one.

Returns:

is_valid (bool) – Whether the date is valid within the current (or newly downloaded) IERS table.

load_data(path=None)[source]

Loads an IERS table in CSV format from path.

update_data(path=None, channel=None)[source]

Update the IERS table, downloading the latest available version.

Time

class coordio.time.Time(value=None, format='jd', scale='UTC')[source]

Bases: object

A time storage and conversion class.

All times, regardless of the input format and scale are converted and stored as TAI Julian dates.

Parameters:
  • value – The input time value. If None, the current datetime.datetime.now system time will be used.

  • format (str) – The format of the input time.

  • scale (str) – The time scale of the input value. Valid values are UTC and TAI.

Variables:

jd (float) – The Julian date in the TAI scale.

get_dut1()[source]

Returns the delta UT1-UTC.

to_now(scale='UTC')[source]

Updates the internal value to the current TAI date.

Parameters:

scale (str) – The scale of the system clock. Assumes it’s UTC but it can be changed if the internal clock is synced to TAI or other timescale.

Returns:

now (float) – The new TAI JD. Also updates the internal value.

to_tdb(longitude=0.0, latitude=0.0, altitude=0.0)[source]

Returns the date converted to JD in the TDB scale.

Parameters:
  • longitude (float) – The East-positive longitude of the site, in degrees.

  • latitude (float) – The site latitude in degrees.

  • altitude (float) – The altitude (elevation) of the site above see level, in meters. Defaults to zero meters.

to_tt()[source]

Returns the date converted to JD in the TT scale.

to_ut1()[source]

Returns the date converted to JD in the UT1 scale.

to_utc()[source]

Returns the date converted to JD in the UTC scale.

property iers

Gets the IERS table.

property jd1

Returns the integer part of the Julian data.

property jd2

Returns the fractional part of the Julian data.

property mjd

Returns the modified Julian date (JD-2400000.5).

coordio.time.utc_to_tai(jd)[source]

Convert UTC to TAI. Wrapper around iauUtctai.

Coordinates

class coordio.coordinate.Coordinate(value, **kwargs)[source]

Bases: ndarray

Base class to represent a coordinate system.

Coordinate is a wrapper around the numpy.ndarray to represent a series of points in a given coordinate system. It handles most of the particularities of subclassing a Numpy array.

A Coordinate instance must be a 2D array of dimensions \(N\times M\) in which N is the number of coordinate points, and M are their axis values. Extra parameters and arrays can be added to the coordinate object. For example, the following code minimally implements an ICRS coordinate system

class ICRS(Coordinate):

__extra_params__ ['name']
__extra_arrays__ = ['epoch', 'pmra', 'pmdec', 'parallax', 'rvel']

Now we can instantiate a set of coordinates

>>> icrs = ICRS([[100., 20.], [200., 10], [300., 5.]],
                pmra=[1, 2, 0], pmdec=[-0.5, 0.5, 0],
                name='my_coordinates')
>>> print(icrs)
ICRS([[100.,  20.],
      [200.,  10.],
      [300.,   5.]])
>>> print(icrs.pmra)
[1, 2, 0]
>>> print(icrs.parallax)
[0, 0, 0]
>>> print(icrs.name)
my_coordinates

Values in __extra_params__ are propagated to the new instance as attributes and can be defined with any user-provided value. __extra_arrays__ on the other hand are expected to be arrays with the same length as the number of coordinates. If the extra array value is not provided when instantiating the coordinate, it is replaced with arrays of zeros of size N (the number of coordinates). If you don’t want the default value to be zeros you’ll need to override the __array_finalize__ method, for example

def __array_finalize__(self, obj):

    super().__array_finalize__(obj)

    if getattr(obj, 'epoch', None) is None:
        self.epoch += 2451545.0

And now the default value for epoch will be 2451545.

Coordinate takes care of slicing the arrays appropriately. If we slice the coordinate along the first axis the resulting object will be another coordinate instance in which both the internal array and the extra ones are sliced correctly

>>> new_icrs = icrs[0:2, :]
>>> print(new_icrs)
ICRS([[100.,  20.],
      [200.,  10.]])
>>> print(icrs.pmra)
[1, 2]

But if we slice it in a way that the resulting array cannot represent the coordinate system anymore, a simple Numpy array will be returned, without any of the additional information

>>> icrs[:, 0]
array([100., 200., 300.])
property coordSysName
class coordio.coordinate.Coordinate2D(value, **kwargs)[source]

Bases: Coordinate

2 dimensional coordinate system

class coordio.coordinate.Coordinate3D(value, **kwargs)[source]

Bases: Coordinate

3 dimensional coordinate system

coordio.coordinate.verifySite(kwargs, strict=True)[source]

Search through kwargs for site parameter, and verify that it looks ok. Otherwise raise a CoordIOError

Parameters:
  • kwargs (dict) – kwargs passed to Coordinate __new__

  • strict (bool) – if True, a time must be specified on the site

coordio.coordinate.verifyWavelength(kwargs, lenArray, strict=True)[source]

Search through kwargs for wavelength parameter. If not existent, return insert GFA wavelength to kwargs (of the right size). If strict and wavelengths do not correspond to Apogee, Boss, or GFA wavelengths, raise a CoordIOError

Parameters:
  • kwargs (dict) – kwargs passed to Coordinate __new__

  • lenArray (int) – length of the array to create if wavelength not present in kwargs

  • strict (bool) – if strict is True, wavelengths may only be those corresponding to Apogee, Boss, or GFA

class coordio.sky.ICRS(value, **kwargs)[source]

Bases: Coordinate2D

A representation of ICRS coordinates.

Parameters:
  • value (numpy.ndarray) – A Nx2 Numpy array with the RA and Dec coordinates of the targets.

  • epoch (numpy.ndarray) – A 1D array with the epoch of the coordinates for each target, as a TDB Julian date (although for most applications the small differences between scales will not matter). Defaults to J2000.

  • pmra (numpy.ndarray) – A 1D array with the proper motion in the RA axis for the N targets, in milliarcsec/yr. Must be a true angle, i.e, it must include the cos(dec) term.

  • pmdec (numpy.ndarray) – A 1D array with the proper motion in the RA axis for the N targets, in milliarcsec/yr.

  • parallax (numpy.ndarray) – A 1D array with the parallax for the N targets, in milliarcsec.

  • rvel (numpy.ndarray) – A 1D array with the radial velocity in km/s, positive when receding.

to_epoch(jd, site=None)[source]

Convert the coordinates to a new epoch.

Parameters:
  • jd (float) – The Julian date, in TAI scale, of the output epoch.

  • site (.Site) – The site of the observation. Used to determine the TDB-TT offset. If not provided, it assumes longitude and latitude zero.

Returns:

icrs (ICRS) – A new ICRS object with the coordinates, proper motion, etc. in the new epoch.

class coordio.sky.Observed(value, **kwargs)[source]

Bases: Coordinate2D

The observed coordinates of a series of targets.

The array contains the Alt/Az coordinates of the targets. Their RA/Dec coordinates can be accessed via the ra and dec attributes. If ICRS or Field is passed, Alt/Az coordinates are computed.

Parameters:
  • value (numpy.ndarray) – A Nx2 Numpy array with the Alt and Az coordinates of the targets, in degrees. Or ICRS instance. Or a Field instance.

  • wavelength (numpy.ndarray) – A 1D array with he observing wavelength, in angstrom. If not explicitly passed, it tries to inheret from value.wavelength, if that doesn’t exist, it is set to default specified in: defaults.wavelength

  • site (.Site) – The site from where observations will occur, along with the time of the observation. Mandatory argument.

Variables:
  • ra (numpy.ndarray) – Nx1 Numpy array, observed RA in degrees

  • dec (numpy.ndarray) – Nx1 Numpy array, observed Dec in degrees

  • ha (numpy.ndarray) – Nx1 Numpy array, hour angle in degrees

  • pa (numpy.ndarray) – Nx1 Numpy array, position angle in degrees. By SOFA: the angle between the direction to the north celestial pole and direction to the zenith. range is [-180, 180]. The sign is according to: -ha –> -pa, +ha –> +pa

classmethod fromEquatorial(coords, hadec=False, site=None, **kwargs)[source]

Initialises Observed coordinates from equatorial topocentric.

Parameters:
  • coords (numpy.ndarray) – Nx2 Numpy array with the RA and Dec coordinates of the targets. These must be topocentric equatorial coordinates that will be converted to horizontal topocentric used to initialise Observed. Alternatively, if hadec=True then the first column in the array must be the hour angle in degrees.

  • hadec (bool) – Whether the coordinates are HA and declination.

  • site (.Site) – The site from where observations will occur, along with the time of the observation. Mandatory argument.

  • kwargs – Other arguments to pass to Observed.

dec: numpy.ndarray
ha: numpy.ndarray
pa: numpy.ndarray
ra: numpy.ndarray
site: Site
wavelength: numpy.ndarray
class coordio.telescope.Field(value, **kwargs)[source]

Bases: Coordinate2D

A representation of Field Coordinates. A spherical coordinate system defined by two angles: theta, phi. Theta is the angle about the optical axis measured from the direction of +RA. Phi is the angle off the optical axis of the telescope. So Phi=0 is the optical axis of the telescope increasing from earth to sky.

In the Cartesian representation on a unit sphere, +x is aligned with +RA, +y is aligned with +Dec, and so +z is aligned with the telescope’s optical axis and points from the earth to sky.

Parameters:
  • value (numpy.ndarray) – A Nx2 array. First column is theta, the azimuthal angle from +RA through +Dec in deg. Second angle is phi, polar angle angle (off axis angle) in deg. Or Observed. Or FocalPlane.

  • field_center (Observed) – An Observed instance containing a single coordinate

Variables:
field_center: Observed
field_warn: numpy.ndarray
x: numpy.ndarray
x_angle: numpy.ndarray
y: numpy.ndarray
y_angle: numpy.ndarray
z: numpy.ndarray
class coordio.telescope.FocalPlane(value, **kwargs)[source]

Bases: Coordinate3D

The focal plane coordinate system is 3D Cartesian with units of mm. The origin is the M1 vertex. +x is aligned with +RA on the image, +y is aligned with +Dec on the image, +z points from earth to sky along the telescope’s optical axis. FocalPlane xy’s are physically rotated 180 degrees with respect to Field xy’s because Cassegrain telescope rotates the image. The choice was made to put +x along +RA on the image, rather than the sky.

For each wavelength, a sphereical surface fit has been determined. The sphere’s center lies along the optical axis at a positive z value. The sphere’s center is considered the ‘ray origin’, so all rays traveling to the focal plane originate from this location.

The ray origin and a wavelength-dependent, z-symmetric, distortion model is used to convert between FocalPlane and Field coordinates. The distortion model is an odd-termed polynomial relating input field angle to output focal plane model.

The ZEMAX models/scripts, and fitting code for this still lives in github.com/csayres/sdssconv, but may be eventually moved to this package.

Parameters:
  • value (numpy.ndarray) – A Nx3 array where [x,y,z] are columns. Or Field. Or Wok.

  • wavelength (float or numpy.ndarray) – A 1D array with he observed wavelength, in angstroms. Currently only values of 5400, 6231, and 16600 are valid, corresponding to BOSS, Apogee, and sdss-r band (GFA). Defaults to GFA wavelength.

  • site (Site) – Used to pick the focal plane model to use which is telescope dependent

  • fpScale (float) – Multiplicative scale factor to apply between wok and focal coords. An adjustment to the focal plane model. Defaults to coordio.defaults.FOCAL_SCALE.

  • use_closest_wavelength (bool) – If True, uses the ZEMAX model with the closes wavelength to the input wavelength.

Variables:
  • b (numpy.ndarray) – A 1D array containing the location of the spherical fit’s center (origin of rays) along the optical axis in mm. Varies with wavelength.

  • R (numpy.ndarray) – 1D array containing the radius of curvature for the spherical fit in mm. Varies with wavelength.

  • field_warn (numpy.ndarray) – boolean array indicating suspect conversions from large field angles

R: numpy.ndarray
b: numpy.ndarray
field_warn: numpy.ndarray
fpScale: float
site: Site
use_closest_wavelength: bool
wavelength: numpy.ndarray
class coordio.wok.Wok(value, **kwargs)[source]

Bases: Coordinate3D

A representation of a Wok Coordinates. A 3D Cartesian coordinate system. The orgin is the vertex of the Wok. -x points toward the Boss slithead. +z points from earth to sky. +y is dictated by right hand rule (so it points toward GFA-S2).

A rotation about focal coords’ +z is applied to observe fields at non-zero position angles

With perfect alignment to the rotator at zero position angle, the xyz axes of the wok are idenentical to those of the focal plane system. The origin of the wok is translated along z with respect to the focal plane system by an amount determined by the ZEMAX models of the focal plane.

True translations and tilts of the wok with respect to the focal plane are measured on sky (at zero position angle) and stored in the etc/wokOrientation.csv file.

Parameters:
  • value (numpy.ndarray) – A Nx3 array where [x,y,z] are columns. Or FocalPlane. Or Tangent.

  • site (Site) – site name determines which wok parameters to use. Mandatory parameter.

  • obsAngle (float) – Position angle of observation. Angle measured from (image) North through East to wok +y. So obsAngle of 45 deg, wok +y points NE. Defaults to zero.

obsAngle: float
site: Site
class coordio.tangent.Tangent(value, **kwargs)[source]

Bases: Coordinate3D

A representation of Tangent Coordinates. A 3D Cartesian coordinate system. The xy plane of this coordinate system is tangent to the wok surface at a specific location (one of holeID in etc/wokCoords.csv). For positioners, +x axis is a aligned with the alpha=0 direction. For GFA’s the -y axis points toward the wok coordinate system origin. +z axis (generally) points from the earth to sky. Origin lies 143 mm above the wok surface. This puts the Tangent xy plane in the plane of the fiber or chip.

Parameters:
  • value (numpy.ndarray) – A Nx3 array where [x,y,z] are columns. Or Wok. Or Positioner. Or Guide.

  • site (Site) – Site name determines which wok parameters to use. Mandatory parameter.

  • holeID (str) – Calid identifier(s), one of calibration.VALID_HOLE_IDS

  • scaleFactor (float) – Multiplicative factor to apply, modeling thermal expansion/contraction of wok holes with respect to each other. Defaults to 1

  • wavelength (float or numpy.ndarray) – wavelength used for projecting rays to tangent surfaces (from sphere model origin). Defaults to GFA wavelength

Variables:
  • xProj (numpy.ndarray) – x projection of coordinate to the xy plane

  • yProj (numpy.ndarray) – y projection of coordinate to xy plane

  • distProj (numpy.ndarray) – Distance of projection (mm) proxy for focus offset positive when above tangent surface, negative when coord is below tangent surface

  • obsAngle (float) – Position angle of observation. Angle measured from (image) North through East to wok +y. So obsAngle of 45 deg, wok +y points NE. Defaults to zero.

distProj: numpy.ndarray
holeID: str
obsAngle: float | numpy.ndarray
scaleFactor: float | numpy.ndarray
site: Site
wavelength: numpy.ndarray
xProj: numpy.ndarray
yProj: numpy.ndarray
class coordio.tangent.TangentNoProj(value, **kwargs)[source]

Bases: Tangent

Class that doesn’t compute projections, intended for internal use only to eliminate a recursion problem.

class coordio.guide.Guide(value, **kwargs)[source]

Bases: Coordinate2D

Guide coordinates are 2D cartesian xy coordinates. Units are pixels. The origin is the lower left corner of the lower left pixel when looking at the CCD through the camera window. Thus guide images are taken in “selfie-mode”. So the coordinate (0.5, 0.5) is the center of the lower left pixel. -y points (generally) toward the wok vertex. True rotations and locations of GFAs in the wok are measured/calibrated after installation and on-sky characterization.

Parameters:
  • value (numpy.ndarray) – A Nx2 array where [x,y] are columns. Or Tangent

  • xBin (int) – CCD x bin factor, defaults to 1

  • yBin (int) – CCD y bin factor, defaults to 1

Variables:

guide_warn (numpy.ndarray) – boolean array indicating coordinates outside the CCD FOV

guide_warn: ndarray
xBin: int
yBin: int
class coordio.guide.GuiderFitter(observatory: str)[source]

Bases: object

Fits guider data, returning the measured offsets.

add_astro(camera: str | int, ra: float, dec: float, obstime: float, x_pixel: float, y_pixel: float)[source]

Adds an astrometric measurement.

add_gimg(path: str | Path, pixels: ndarray | None = None)[source]

Processes a raw GFA image, runs astrometry.net, adds the solution.

add_proc_gimg(path: str | Path, pixels: ndarray | None = None)[source]

Adds an astrometric solution from a proc-gimg image.

add_wcs(camera: str | int, wcs: WCS, obstime: float, pixels: ndarray | None = None)[source]

Adds a camera measurement from a WCS solution.

calculate_rms(gfa_wok: DataFrame, astro_wok: DataFrame, scale: float = 1, cameras: list[int] | None = None)[source]

Calculates the RMS of the measurements.

Parameters:
  • gfa_wok – GFA to wok coordinates data. The data frame contains three columns: gfa_id, xwok, and ywok.

  • astro_wok – Astrometric solution to wok coordinates data. The data frame contains three columns: gfa_id, xwok, and ywok.

  • scale – Factor by which to scale coordinates.

  • cameras – The cameras to use to calculate the global RMS.

Returns:

rms – A data frame with columns gfa_id, xrms, yrms, rms for the x, y, and combined RMS measurements for each camera. An additional gfa_id=0 is added with the RMS measurements for all GFA cameras combined. RMS values are returned in mm on wok coordinates.

fit(field_ra: float, field_dec: float, field_pa: float = 0.0, offset: tuple[float, float, float] | ndarray = (0.0, 0.0, 0.0), scale_rms: bool = False, only_radec: bool = False, cameras: list[int] | None = None)[source]

Performs the fit and returns a GuiderFit object.

The fit is performed using umeyama.

Parameters:
  • field_ra – The field centre RA.

  • field_dec – The field centre declination.

  • field_pa – The field centre position angle.

  • offset – The offset in RA, Dec, PA to apply.

  • scale_rms – Whether to correct the RMS using the measured scale factor.

  • only_radec – If True, fits only translation. Useful when only one or two cameras are solving. In this case the rotation and scale offsets will still be calculated using the Umeyama model, but the ra and dec offsets are overridden with a simple translation. The GuiderFit object will have only_radec=True.

plot_fit(filename: str | Path)[source]

Plot the fit results.

reset()[source]

Clears the astrometric data.

coordio.guide.cross_match(measured_xy: ndarray, reference_xy: ndarray, reference_radec: ndarray, x_size: int, y_size: int, cross_corrlation_shift: bool = True, blur: float = 5, distance_upper_bound: int = 10, **kwargs)[source]

Determines the shift between two sets of points using cross-correlation.

Constructs 2D images from reference and measured data sets and calculates the image shift using cross-correlation registration. It then associated reference to measured detections using nearest neighbours and builds a WCS using the reference on-sky positions.

Parameters:
  • measured_xy – A 2D array with the x and y coordinates of each measured point.

  • reference_xy – A 2D array with the x and y coordinates of each reference point.

  • reference_radec – A 2D array with the ra and dec coordinates of each reference point.

  • x_size – Size of the image for 2D cross-correlation.

  • y_size – Size of the image for 2D cross-correlation.

  • blur – The sigma, in pixels, of the Gaussian kernel used to convolve the images.

  • distance_upper_bound – Maximum distance, in pixels, for KD tree nearest neighbours.

  • kwargs – Other arguments to pass to skimage.registration.phase_cross_correlation.

Returns:

wcs – A tuple with the WCS of the solution and the translation invariant normalized RMS error between reference and moving image (see phase_cross_correlation).

coordio.guide.gfa_to_radec(observatory: str, x_pix: float, y_pix: float, gfa_id: int, bore_ra: float, bore_dec: float, position_angle: float = 0, offset_ra: float = 0, offset_dec: float = 0, offset_pa: float = 0, obstime: float | None = None, icrs: bool = False)[source]

Converts from a GFA pixel to observed RA/Dec. Offsets are in arcsec.

coordio.guide.radec_to_gfa(observatory: str, ra: float | ndarray, dec: float | ndarray, gfa_id: int, bore_ra: float, bore_dec: float, position_angle: float = 0, offset_ra: float = 0, offset_dec: float = 0, offset_pa: float = 0, obstime: float | None = None)[source]

Converts from observed RA/Dec to GFA pixel, taking into account offsets.

coordio.guide.umeyama(X, Y)[source]

Rigid alignment of two sets of points in k-dimensional Euclidean space.

Given two sets of points in correspondence, this function computes the scaling, rotation, and translation that define the transform TR that minimizes the sum of squared errors between TR(X) and its corresponding points in Y. This routine takes O(n k^3)-time.

Parameters:
  • X – A k x n matrix whose columns are points.

  • Y – A k x n matrix whose columns are points that correspond to the points in X

Returns:

  • c,R,t – The scaling, rotation matrix, and translation vector defining the linear map TR as

    TR(x) = c * R * x + t
    

    such that the average norm of TR(X(:, i) - Y(:, i)) is minimized.

  • Copyright (Carlo Nicolini, 2013)

  • Code adapted from the Mark Paskin Matlab version from

  • http (//openslam.informatik.uni-freiburg.de/data/svn/tjtf/trunk/matlab/ralign.m)

  • See paper by Umeyama (1991)

class coordio.positioner.PositionerApogee(value, **kwargs)[source]

Bases: PositionerBase

class coordio.positioner.PositionerBase(value, **kwargs)[source]

Bases: Coordinate2D

A representation of Positioner coordinates. Alpha/Beta coordinates in degrees.

When converting from tangent coordinates, Robot’s are always “right handed” with alpha in [0, 360], beta in [0, 180]

Parameters:
  • value (numpy.ndarray) – A Nx2 array with columns [alpha, beta] in degrees

  • site (Site) – mandatory parameter

  • holeID (str) – valid identifier, one of calibration.VALID_HOLE_IDS

holeID: str
positioner_warn: numpy.ndarray
site: Site
class coordio.positioner.PositionerBoss(value, **kwargs)[source]

Bases: PositionerBase

class coordio.positioner.PositionerMetrology(value, **kwargs)[source]

Bases: PositionerBase

Conversions

coordio.conv.cart2FieldAngle(x, y, z)[source]

Convert cartesian point on unit sphere of sky to (ZEMAX-style) field angles in degrees.

Zemax defines field angles like this: Positive field angles imply positive slope for the ray in that direction, and thus refer to negative coordinates on distant objects. ZEMAX converts x field angles ( αx ) and y field angles ( αy ) to ray direction cosines using the following formulas:

tanαx = l/n tanαy = m/n l^2 + m^2 + n^2 = 1

where l, m, and n are the x, y, and z direction cosines.

Parameters:
  • x (scalar or 1D array)

  • y (scalar or 1D array)

  • z (scalar or 1D array)

Returns:

result (list) – [xField, yField] ZEMAX style field coordinates (degrees)

coordio.conv.cart2Sph(x, y, z)[source]

Convert cartesian coordinates to spherical coordinates theta, phi in degrees

phi is polar angle measure from z axis theta is azimuthal angle measured from x axis

Parameters:
  • x (scalar or 1D array)

  • y (scalar or 1D array)

  • z (scalar or 1D array)

Returns:

result (list) – [theta, phi] degrees

coordio.conv.fieldAngle2Cart(xField, yField)[source]

Convert (ZEMAX-style) field angles in degrees to a cartesian point on the unit sphere

Zemax defines field angles like this: Positive field angles imply positive slope for the ray in that direction, and thus refer to negative coordinates on distant objects. ZEMAX converts x field angles ( αx ) and y field angles ( αy ) to ray direction cosines using the following formulas:

tanαx = l/n tanαy = m/n l^2 + m^2 + n^2 = 1

where l, m, and n are the x, y, and z direction cosines.

Parameters:
  • xField (scalar or 1D array) – zemax style x field

  • yField (scalar or 1D array) – zemax style y field

Returns:

result (list) – [x,y,z] coordinates on unit sphere

coordio.conv.fieldToFocal(thetaField, phiField, site, waveCat)[source]

Convert spherical field coordinates to the modeled spherical position on the focal plane. Origin is the M1 vertex.

Raises warnings for coordinates at large field angles.

Parameters:
  • thetaField (float or numpy.ndarray) – azimuthal angle of field coordinate (deg)

  • phiField (float or numpy.ndarray) – polar angle of field coordinate (deg)

  • site (str) – “APO” or “LCO”

  • waveCat (str) – “Boss”, “Apogee”, or “GFA”

  • Result

  • -------

  • xFocal (float or numpy.ndarray) – spherical x focal coord mm

  • yFocal (float or numpy.ndarray) – spherical y focal coord mm

  • zFocal (float or numpy.ndarray) – spherical z focal coord mm

  • R (float) – radius of curvature for sphere

  • b (float) – z location of center of sphere

  • fieldWarn (bool or boolean array) – True if angle off-axis is large enough to be suspicious

coordio.conv.fieldToObserved(x, y, z, altCenter, azCenter, pa)[source]

Convert xyz unit-spherical field coordinates to Alt/Az

the inverse of observetToField

Az=0 is north, Az=90 is east

+x points along +RA, +y points along +Dec, +z points along boresight from telescope toward sky

Parameters:
  • x (scalar or 1D array) – unit-spherical x field coord (aligned with +RA)

  • y (scalar or 1D array) – unit-spherical y field coord (aligned with +Dec)

  • z (scalar or 1D array) – unit-spherical z coord

  • altCenter (float) – altitude in degrees, field center

  • azCenter (float) – azimuth in degrees, field center

  • pa (float) – position angle of field center

Returns:

  • alt (float or 1D array) – altitude of coordinate in degrees

  • az (float or 1D array)

coordio.conv.focalToField(xFocal, yFocal, zFocal, site, waveCat)[source]

Convert xyz focal position to unit-spherical field position. xyzFocal do not need to lie on the modeled sphere, and the spherical radius of curvature is not needed, only the origin of the sphere is needed (b).

Raises warnings for coordinates at large field angles.

Parameters:
  • xFocal (float or 1D array) – x focal coord mm

  • yFocal (float or 1D array) – y focal coord mm

  • zFocal (float or 1D array) – z focal coord mm, +z points along boresight toward sky.

  • site (str) – “APO” or “LCO”

  • waveCat (str) – “Boss”, “Apogee”, or “GFA”

  • Result

  • -------

  • thetaField (float or 1D array) – azimuthal field coordinate (deg)

  • phiField (float or 1D array) – polar field coordinate (deg)

  • fieldWarn (bool or boolean array) – True if angle off-axis is large enough to be suspicious

coordio.conv.focalToWok(xFocal, yFocal, zFocal, positionAngle=0, xOffset=0, yOffset=0, zOffset=0, tiltX=0, tiltY=0, fpScale=1, projectFlat=True)[source]

Convert xyz focal coordinates in mm to xyz wok coordinates in mm.

The origin of the focal coordinate system is the M1 vertex. focal +y points toward North, +x points toward E. The origin of the wok coordinate system is the wok vertex. -x points toward the boss slithead. +z points from telescope to sky.

Tilt is applied about x axis then y axis.

Note: currently zFocal irrelevant using a flat wok model.

Parameters:
  • xFocal (scalar or 1D array) – x position of object on focal plane mm (+x aligned with +RA on image)

  • yFocal (scalar or 1D array) – y position of object on focal plane mm (+y aligned with +Dec on image)

  • zFocal (scalar or 1D array) – z position of object on focal plane mm (+z aligned boresight and increases from the telescope to the sky)

  • positionAngle (scalar) – position angle deg. Angle measured from (image) North through East to wok +y. So position angle of 45 deg, wok +y points NE

  • site (str) – either “APO” or “LCO”

  • xOffset (scalar or None) – x position (mm) of wok origin (vertex) in focal coords calibrated

  • yOffset (scalar) – y position (mm) of wok origin (vertex) in focal coords calibratated

  • zOffset (scalar) – z position (mm) of wok origin (vertex) in focal coords calibratated

  • tiltX (scalar) – tilt (deg) of wok about focal x axis at PA=0 calibrated

  • tiltY (scalar) – tilt (deg) of wok about focal y axis at PA=0 calibrated

  • fpScale (scalar) – radial scale multiplier

  • projectFlat (bool) – if true, assume flat wok model

Returns:

  • xWok (scalar or 1D array) – x position of object in wok space mm (+x aligned with +RA on image)

  • yWok (scalar or 1D array) – y position of object in wok space mm (+y aligned with +Dec on image)

  • zWok (scalar or 1D array) – z position of object in wok space mm (+z aligned boresight and increases from the telescope to the sky)

coordio.conv.guideToTangent(xPix, yPix, xBin=1, yBin=1, x_0=0, x_1=0, y_0=0, y_1=0)[source]
coordio.conv.observedToField(alt, az, altCenter, azCenter, pa)[source]

Convert Az/Alt coordinates to cartesian coords on the unit sphere with the z axis aligned with field center (boresight)

Az=0 is north, Az=90 is east

Resulting Cartesian coordinates have the following convention +x points along +RA, +y points along +Dec

Parameters:
  • alt (float or 1D array) – altitude in degrees, positive above horizon, negative below

  • az (float or 1D array) – azimuth in degrees az=0 is north az=90 is east

  • altCenter (float) – altitude in degrees, field center

  • azCenter (float) – azimuth in degrees, field center

  • pa (float) – position angle of field center

Returns:

  • theta (float or 1D array) – azimuthal angle of field coordinate (deg)

  • phi (float or 1D array) – polar angle of field coordinat (deg off axis)

coordio.conv.positionerToTangent(alphaDeg, betaDeg, xBeta, yBeta, la=7.4, alphaOffDeg=0, betaOffDeg=0)[source]

Determine tangent coordinates (mm) of xBeta, yBeta coords in mm from alpha/beta angle. CPP version.

todo: include hooks for positioner non-linearity

Parameters:
  • alphaDeg (scalar or 1D array) – alpha angle in degrees

  • betaDeg (scalar or 1D array) – beta angle in degrees

  • xBeta (scalar or 1D array) – x position (mm) in beta arm frame

  • yBeta (scalar or 1D array) – y position (mm) in beta arm frame

  • la (scalar or 1D array) – length (mm) of alpha arm

  • alphaOffDeg (scalar) – alpha zeropoint offset (deg)

  • betaOffDeg (scalar) – beta zeropoint offset (deg)

Returns:

  • xTangent (scalar or 1D array) – x position (mm) in tangent coordinates

  • yTangent (scalar or 1D array) – y position (mm) in tangent coordinates

coordio.conv.positionerToWok(alphaDeg, betaDeg, xBeta, yBeta, la, alphaOffDeg, betaOffDeg, basePos, iHat, jHat, kHat, dx, dy, elementHeight=143.1)[source]

alpha/betaDeg: angular coords for arm angles (deg) xyBeta: location of fiber in beta arm (mm) la: alpha arm length alpha/betaOffDeg: calibrated angular zeropoint offsets for arm b: xyz locations for holes in the wok ijkHat: unit vectors to capture any calibrated rotation dxy: offsets of element from hole centered (calibrated) elementHeight (vertical distance from wok surface to fiber surface)

returns xWok, yWok, zWok

coordio.conv.proj2XYplane(x, y, z, rayOrigin)[source]

Given a point x, y, z orginating from rayOrigin, project this point onto the xy plane.

Parameters:
  • x (scalar or 1D array) – x position of point to be projected

  • y (scalar or 1D array) – y position of point to be projected

  • z (scalar or 1D array) – z position of point to be projected

  • rayOrigin (3-vector) – [x,y,z] origin of ray(s)

Returns:

  • xProj (scalar or 1D array) – projected x value

  • yProj (scalar or 1D array) – projected y value

  • zProj (scalar or 1D array) – projected z value (should be zero always!)

  • projDist (scalar or 1D array) – projected distance (a proxy for something like focus offset)

coordio.conv.repeat_scalar(input, n)[source]

Returns an array with an scalar repeated n times. Does not affect arrays.

coordio.conv.sph2Cart(theta, phi, r=1)[source]

Convert spherical coordinates theta, phi in degrees to cartesian coordinates on unit sphere.

phi is polar angle measure from z axis theta is azimuthal angle measured from x axis

Parameters:
  • theta (scalar or 1D array) – degrees, azimuthal angle

  • phi (scalar or 1D array) – degrees, polar angle

  • r (scalar) – radius of curvature. Default to 1 for unit sphere

Returns:

result (list) – [x,y,z] coordinates on sphere

coordio.conv.tangentToGuide(xTangent, yTangent, xBin=1, yBin=1, x_0=0, x_1=0, y_0=0, y_1=0)[source]
coordio.conv.tangentToPositioner(xTangent, yTangent, xBeta, yBeta, la=7.4, alphaOffDeg=0, betaOffDeg=0, lefthand=False)[source]

Determine alpha/beta positioner angles that place xBeta, yBeta coords in mm at xTangent, yTangent. CPP Version.

todo: include hooks for positioner non-linearity

Parameters:
  • xTangent (scalar or 1D array) – x position (mm) in tangent coordinates

  • yTangent (scalar or 1D array) – y position (mm) in tangent coordinates

  • xBeta (scalar or 1D array) – x position (mm) in beta arm frame

  • yBeta (scalar or 1D array) – y position (mm) in beta arm frame

  • la (scalar or 1D array) – length (mm) of alpha arm

  • alphaOffDeg (scalar) – alpha zeropoint offset (deg)

  • betaOffDeg (scalar) – beta zeropoint offset (deg)

  • lefthand (boolean) – whether to return a lefthand or righthand parity robot

Returns:

  • alphaDeg (scalar or 1D array) – alpha angle in degrees

  • betaDeg (scalar or 1D array) – beta angle in degrees

  • isOK (boolean or 1D boolean array) – True if point physically accessible, False otherwise

coordio.conv.tangentToWok(xTangent, yTangent, zTangent, b, iHat, jHat, kHat, elementHeight=143.1, scaleFac=1, dx=0, dy=0, dz=0)[source]

Convert from tangent coordinates at b to wok coordinates.

xyz Wok coords are mm with origin at wok vertex. +z points from wok toward M2. -x points toward the boss slithead

In the tangent coordinate frame the xy plane is tangent to the wok surface at xyz position b. The origin is set to be elementHeight above the wok surface.

Parameters:
  • xTangent (scalar or 1D array) – x position (mm) in tangent coordinates

  • yTangent (scalar or 1D array) – y position (mm) in tangent coordinates

  • zTangent (scalar or 1D array) – z position (mm) in tangent coordinates

  • b (3-vector) – x,y,z position (mm) of element hole on wok surface measured in wok coords

  • iHat (3-vector) – x,y,z unit vector in wok coords that indicate the direction of the tangent coordinate x axis

  • jHat (3-vector) – x,y,z unit vector in wok coords that indicate the direction of the tangent coordinate y axis

  • kHat (3-vector) – x,y,z unit vector in wok coords that indicate the direction of the tangent coordinate z axis

  • elementHeight (scalar) – height (mm) of positioner/GFA chip above wok surface

  • scaleFac (scalar) – scale factor to apply to b to account for thermal expansion of wok. scale is applied to b radially

  • dx (scalar) – x offset (mm), calibration to capture small displacements of tangent x

  • dy (scalar) – y offset (mm), calibration to capture small displacements of tangent y

  • dz (scalar) – z offset (mm), calibration to capture small displacements of tangent x

Returns:

  • xWok (scalar or 1D array) – x position of object in wok space mm (+x aligned with +RA on image at PA = 0)

  • yWok (scalar or 1D array) – y position of object in wok space mm (+y aligned with +Dec on image at PA = 0)

  • zWok (scalar or 1D array) – z position of object in wok space mm (+z aligned boresight and increases from the telescope to the sky)

coordio.conv.wok2mirror(x, y, ipa)[source]

Convert to “mirror coordinates” in mm

coordio.conv.wokToFocal(xWok, yWok, zWok, positionAngle=0, xOffset=0, yOffset=0, zOffset=0, tiltX=0, tiltY=0, fpScale=1, b=None, R=None)[source]

Convert xyz wok coordinates in mm to xyz focal coordinates in mm.

The origin of the focal coordinate system is the M1 vertex. focal +y points toward North, +x points toward E. The origin of the wok coordinate system is the wok vertex. -x points toward the boss slithead. +z points from telescope to sky.

Tilt is applied first about x axis then y axis.

Note: currently zWok is irrelevant using a flat wok model.

Parameters:
  • xWok (scalar or 1D array) – x position of object in wok space mm (+x aligned with +RA on image at PA=0)

  • yWok (scalar or 1D array) – y position of object in wok space mm (+y aligned with +Dec on image at PA=0)

  • zWok (scalar or 1D array) – z position of object in wok space mm (+z aligned boresight and increases from the telescope to the sky)

  • positionAngle (scalar) – position angle deg. Angle measured from (image) North through East to wok +y. So position angle of 45 deg, wok +y points NE

  • xOffset (scalar or None) – x position (mm) of wok origin (vertex) in focal coords calibrated

  • yOffset (scalar) – y position (mm) of wok origin (vertex) in focal coords calibratated

  • zOffset (scalar) – z position (mm) of wok origin (vertex) in focal coords calibratated

  • tiltX (scalar) – tilt (deg) of wok about focal x axis at PA=0 calibrated

  • tiltY (scalar) – tilt (deg) of wok about focal y axis at PA=0 calibrated

  • fpScale (scalar) – radial scale multiplier

  • b (float) – z location of the center of curvature, dist from M1 vertex (mm)

  • R (float) – radius of curvature (mm)

Returns:

  • xFocal (scalar or 1D array) – x position of object in focal coord sys mm (+x aligned with +RA on image)

  • yFocal (scalar or 1D array) – y position of object in focal coord sys mm (+y aligned with +Dec on image)

  • zFocal (scalar or 1D array) – z position of object in focal coord sys mm (+z aligned boresight and increases from the telescope to the sky)

coordio.conv.wokToPositioner(xWok, yWok, zWok, xBeta, yBeta, la, alphaOffDeg, betaOffDeg, basePos, iHat, jHat, kHat, dx, dy, elementHeight=143.1)[source]

xyzWok: xyz wok location of fiber xyBeta: location of fiber in beta arm (mm) la: alpha arm length alpha/betaOffDeg: calibrated angular zeropoint offsets for arm b: xyz locations for holes in the wok ijkHat: unit vectors to capture any calibrated rotation dxy: offsets of element from hole centered (calibrated) elementHeight (vertical distance from wok surface to fiber surface)

returns alpha, beta angles (deg)

coordio.conv.wokToTangent(xWok, yWok, zWok, b, iHat, jHat, kHat, elementHeight=143.1, scaleFac=1, dx=0, dy=0, dz=0)[source]

Convert an array of wok coordinates to tangent coordinates.

xyz Wok coords are mm with orgin at wok vertex. +z points from wok toward M2. -x points toward the boss slithead

In the tangent coordinate frame the xy plane is tangent to the wok surface at xyz position b. The origin is set to be elementHeight above the wok surface.

Parameters:
  • xWok (1D array) – x position of object in wok space mm (+x aligned with +RA on image at PA = 0)

  • yWok (1D array) – y position of object in wok space mm (+y aligned with +Dec on image at PA = 0)

  • zWok (1D array) – z position of object in wok space mm (+z aligned boresight and increases from the telescope to the sky)

  • b (3-element 1D vector) – x,y,z position (mm) of each hole element on wok surface measured in wok coords

  • iHat (3-element 1D vector) – x,y,z unit vector in wok coords that indicate the direction of the tangent coordinate x axis for each hole.

  • jHat (3-element 1D vector) – x,y,z unit vector in wok coords that indicate the direction of the tangent coordinate y axis for each hole.

  • kHat (3-element 1D vector) – x,y,z unit vector in wok coords that indicate the direction of the tangent coordinate z axis for each hole.

  • elementHeight (1D array) – height (mm) of positioner/GFA chip above wok surface

  • scaleFac (1D array) – scale factor to apply to b to account for thermal expansion of wok. scale is applied to b radially

  • dx (1D array) – x offset (mm), calibration to capture small displacements of tangent x

  • dy (1D array) – y offset (mm), calibration to capture small displacements of tangent y

  • dz (1D array) – z offset (mm), calibration to capture small displacements of tangent x

Returns:

  • xTangent (1D array) – x position (mm) in tangent coordinates

  • yTangent (1D array) – y position (mm) in tangent coordinates

  • zTangent (1D array) – z position (mm) in tangent coordinates