auromat.coordinates.transform module

This module contains fast and memory-efficient algorithms to convert coordinates between different reference frames or representations.

Note that some functions depend on the IGRF model whose parameters are defined in the igrf module. In this version, parameters are defined until the year 2015.

auromat.coordinates.transform.spherical_to_cartesian(r, lat, lon, astuple=True)[source]

Convert spherical to cartesian coordinates. Inputs must be arrays.

Equivalent to astropy.coordinates.distances.spherical_to_cartesian for array inputs but uses less memory and is faster.

Return type:tuple (x,y,z) of ndarray’s with shape as input
auromat.coordinates.transform.cartesian_to_spherical(x, y, z, with_radius=True)[source]

Convert cartesian to spherical coordinates. Inputs must be arrays.

Equivalent to astropy.coordinates.distances.cartesian_to_spherical for array inputs but uses less memory and is faster.

Return type:tuple (r,lat,lon) or (lat,lon) of ndarray’s with shape as input
auromat.coordinates.transform.geodetic2Ecef(lat, lon, h, a=6378.137, b=6356.752314245179)[source]

Converts geodetic to Earth Centered, Earth Fixed coordinates.

Parameters h, a, and b must be given in the same unit. The values of the return tuple then also have this unit.

Parameters:
  • lat – latitude(s) in radians
  • lon – longitude(s) in radians
  • h – height(s)
  • a – equatorial axis of the ellipsoid of revolution
  • b – polar axis of the ellipsoid of revolution
Return type:

tuple (x,y,z)

auromat.coordinates.transform.geodetic2EcefZero(lat, lon, a=6378.137, b=6356.752314245179)[source]

Fast version of geodetic2Ecef() for h=0.

Parameters:
  • lat – latitude(s) in radians
  • lon – longitude(s) in radians
  • a – equatorial axis of the ellipsoid of revolution
  • b – polar axis of the ellipsoid of revolution
Return type:

tuple (x,y,z)

auromat.coordinates.transform.ecef2Geodetic(x, y, z, a=6378.137, b=6356.752314245179)[source]

Convert ECEF to geodetic coordinates.

This function uses the Bowring algorithm from 1985.

The accuracy is at least 11 decimals (in degrees).

Parameters:
  • x,y,z – ECEF coordinates
  • a – equatorial axis of the ellipsoid of revolution
  • b – polar axis of the ellipsoid of revolution
Return type:

tuple (lat,lon) in radians

auromat.coordinates.transform.j2000ToLatLon(j2000Vecs, time_)[source]

Convert cartesian J2000 coordinates to geodetic coordinates.

Parameters:
  • j2000Vecs – shape (n,3)
  • time (datetime) –
Return type:

tuple (latitudes, longitudes) in degrees

auromat.coordinates.transform.latLonToJ2000(lat, lon, h, time_)[source]

Convert geodetic coordinates to cartesian J2000 coordinates.

Parameters:
  • lat – scalar or 1d-array, degrees
  • lon – scalar or 1d-array, degrees
  • h – scalar or 1d-array
  • time (datetime) –
Return type:

tuple (latitudes, longitudes) in degrees

auromat.coordinates.transform.smLonToMLT(smlons, out=None)[source]

Convert solar magnetic longitudes to magnetic local time.

Parameters:smlons – in degrees [-180,180]
Returns:magnetic local time [0,24]
auromat.coordinates.transform.mltToSmLon(mlt, out=None)[source]

Convert magnetic local time to solar magnetic longitudes.

Parameters:mlt – in hours [0,24]
Returns:solar magnetic longitudes in degrees [-180,180]
auromat.coordinates.transform.j2000ToMLatMLT(j2000Vecs, time_)[source]

Convert cartesian J2000 coordinates to MLat/MLT coordinates using the IGRF model.

MLat = geomagnetic latitude

MLT = magnetic local time

Parameters:
  • j2000Vecs – shape (n,3)
  • time (datetime) –
Return type:

tuple (mlat, mlt) in (degrees,hours)

auromat.coordinates.transform.geoToMLatMLT(geoVecs, time_)[source]

Convert ECEF coordinates to MLat/MLT coordinates using the IGRF model.

MLat = geomagnetic latitude

MLT = magnetic local time

Parameters:
  • geoVecs – shape (n,3)
  • time (datetime) –
Return type:

tuple (mlat, mlt) in (degrees,hours)

auromat.coordinates.transform.smToLatLon(smlats, smlons, time_)[source]

Convert solar magnetic to geodetic coordinates using the IGRF model.

Parameters:
  • smlats – in degrees [-90,90]
  • smlons – in degrees [-180,180]
  • time
Return type:

tuple (latitudes, longitudes) in degrees

auromat.coordinates.transform.northGeomagneticPoleLocation(date)[source]

Calculate the approximate position of the north geomagnetic pole for the given date using the IGRF model.

Parameters:date (datetime) –
Return type:named (latitude, longitude) tuple, in degrees