Orthogonal Magnetic Coordinate Systems

coord_transforms.cd2geodetic(lat_cd_arr, lon_cd_arr, alt_cd_arr, decimals=3, year=2021.0)

Transformation from Centered Dipole (CD) (lat, lon, alt) to Geodetic (lat, lon, alt).

Author: Giorgio Savastano (giorgiosavastano@gmail.com)

Parameters
  • lat_cd_arr (np.ndarray) – array containing latitude component of CD coordinates in degrees

  • lon_cd_arr (np.ndarray) – array containing longitude component of CD coordinates in degrees

  • alt_cd_arr (np.ndarray) – array containing altitude (N.B. not the radial distance) component of CD coordinates in meters

  • decimals (int, default=3) – Number of decimal places to round to. If decimals is negative, it specifies the number of positions to the left of the decimal point

  • year (float, default=2021.0) – year for computing the IGRF Gauss coefficients

Return type

tuple[np.ndarray, np.ndarray, np.ndarray]

coord_transforms.compute_pos_centered_dipole_north_pole(year=2021.0)

Function that returns the position of the centered dipole (CD) northern pole using first three Gauss coefficients. Coordinates of the CD northern pole are returned in Geocentric Spherical (colat, long) coordinates [rad].

Author: Giorgio Savastano (giorgiosavastano@gmail.com)

Parameters

year (float, default=2021.0) – year for the computation

Return type

tuple[float, float]

coord_transforms.eccdf2ecef(x_cd, y_cd, z_cd, year=2021.0)

Transformation from Centered Dipole (CD) Cartesian (x, y, z) to Geocentric Cartesian (x, y, z) coordinates.

Author: Giorgio Savastano (giorgiosavastano@gmail.com)

Parameters
  • x_cd (np.ndarray) – array containing x component of CD coordinates

  • y_cd (np.ndarray) – array containing y component of CD coordinates

  • z_cd (np.ndarray) – array containing z component of CD coordinates

  • year (float, default=2021.0) – year for computing the IGRF Gaus coefficients

Return type

tuple[np.ndarray, np.ndarray, np.ndarray]

coord_transforms.ecef2eccdf(x_geoc, y_geoc, z_geoc, year=2021.0)

Transformation from Geocentric Cartesian (x, y, z) to Centered Dipole (CD) Cartesian (x, y, z) coordinates.

Author: Giorgio Savastano (giorgiosavastano@gmail.com)

Parameters
  • x_geoc (np.ndarray) – array containing x coordinates

  • y_geoc (np.ndarray) – array containing y coordinates

  • z_geoc (np.ndarray) – array containing z coordinates

  • year (float, default=2021.0) – year for computing the IGRF Gaus coefficients

Returns

tuple

Return type

tuple[np.ndarray, np.ndarray, np.ndarray]

coord_transforms.ecef2spherical(x, y, z, decimals=2)

Transformation from geocentric Cartesian (x, y, z) to Spherical (colat, long, r) coordinates.

Author: Giorgio Savastano (giorgiosavastano@gmail.com)

Parameters
  • x (np.ndarray) – array containing x coordinates

  • y (np.ndarray) – array containing y coordinates

  • z (np.ndarray) – array containing z coordinates

  • decimals (int, default=2) – Number of decimal places to round to. If decimals is negative, it specifies the number of positions to the left of the decimal point

Return type

tuple[np.ndarray, np.ndarray, np.ndarray]

coord_transforms.geodetic2cd(gglat_deg_array, gglon_deg_array, ggalt_km_array, decimals=2, year=2021.0)

Transformation from Geodetic (lat, lon, alt) to Centered Dipole (CD) (lat, lon, alt).

Author: Giorgio Savastano (giorgiosavastano@gmail.com)

Parameters
  • gglon_deg_array (np.ndarray) – array containing geodetic longitude values in degrees

  • gglat_deg_array (np.ndarray) – array containing geodetic latitude values in degrees

  • ggalt_km_array (np.ndarray) – array containing geodetic altitude values in km

  • decimals (int, default=2) – Number of decimal places to round to. If decimals is negative, it specifies the number of positions to the left of the decimal point.

  • year (float, default=2021.0) – year for computing the IGRF Gauss coefficients

Returns

CD lat, lon, alt arrays

Return type

tuple[np.ndarray, np.ndarray, np.ndarray]

coord_transforms.geodetic2ecef(geod_lati_deg, geod_long_deg, geod_alti_m)

Conversion from Geodetic (lat, lon, alt) to geocentric Cartesian (x, y, z) coordinates.

Author: Giorgio Savastano (giorgiosavastano@gmail.com)

Parameters
  • geod_long_deg (np.ndarray) – array containing geodetic longitudes in degrees

  • geod_lati_deg (np.ndarray) – array containing geodetic longitudes in degrees

  • geod_alti_m (np.ndarray) – array containing geodetic altitude in meters

Return type

tuple[np.ndarray, np.ndarray, np.ndarray]

coord_transforms.mlon2mlt(mlon_cd, datetime, ssheight=318550)

Computes the magnetic local time at the specified magnetic longitude and UT.

Author: Giorgio Savastano (giorgiosavastano@gmail.com) adapted from apexpy.

Parameters
  • mlon_cd (array_like) – Centered Dipole (CD) magnetic longitude

  • datetime (datetime.datetime) – Date and time

  • ssheight (float, optional) – Altitude in km to use for converting the subsolar point from geographic to magnetic coordinates. A high altitude is used to ensure the subsolar point is mapped to high latitudes, which prevents the South-Atlantic Anomaly (SAA) from influencing the MLT.

Returns

Magnetic local time [0, 24)

Return type

ndarray or float

Notes

To compute the MLT, we find the apex longitude of the subsolar point at the given time. Then the MLT of the given point will be computed from the separation in magnetic longitude from this point (1 hour = 15 degrees).

coord_transforms.spherical2ecef(colat_geoc_arr, long_geoc_arr, radial_dist_geoc_arr)

Transformation from Geocentric (GEO) Spherical (colat, lon, r) to Cartesian (x,y,z).

Author: Giorgio Savastano (giorgiosavastano@gmail.com)

Parameters
  • colat_geoc_arr (np.ndarray) – array containing colatitude component of Geocentric coordinates in rad

  • long_geoc_arr (np.ndarray) – array containing longitude component of Geocentric coordinates in rad

  • radial_dist_geoc_arr (np.ndarray) – array containing radial distance component of Geocentric coordinates in m