Calculates coordinate conversions, field-line mapping, and base vectors.
Parameters:
date (NoneType, float, dt.date, or dt.datetime, optional) – Determines which IGRF coefficients are used in conversions. Uses
current date as default. If float, use decimal year. If None, uses
current UTC. (default=None)
refh (float, optional) – Reference height in km for apex coordinates, the field lines are mapped
to this height. (default=0)
datafile (str or NoneType, optional) – Path to custom coefficient file, if None uses apexsh.dat file
(default=None)
fortranlib (str or NoneType, optional) – Path to Fortran Apex CPython library, if None uses linked library file
(default=None)
Variables:
year (float) – Decimal year used for the IGRF model
RE (float) – Earth radius in km, defaults to mean Earth radius
refh (float) – Reference height in km for apex coordinates
datafile (str) – Path to coefficient file
fortranlib (str) – Path to Fortran Apex CPython library
igrf_fn (str) – IGRF coefficient filename
Notes
The calculations use IGRF-13 with coefficients from 1900 to 2025 [1].
dest (str) – Output coordinate system, accepts ‘geo’, ‘apex’, ‘qd’, or ‘mlt’
height (array_like, optional) – Altitude in km
datetime (datetime.datetime) – Date and time for MLT conversions (required for MLT conversions)
precision (float, optional) – Precision of output (degrees) when converting to geo. A negative
value of this argument produces a low-precision calculation of
geodetic lat/lon based only on their spherical harmonic
representation.
A positive value causes the underlying Fortran routine to iterate
until feeding the output geo lat/lon into geo2qd (APXG2Q) reproduces
the input QD lat/lon to within the specified precision (all
coordinates being converted to geo are converted to QD first and
passed through APXG2Q).
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:
lat (ndarray or float) – Converted latitude (if converting to MLT, output latitude is apex)
lon (ndarray or float) – Converted longitude or MLT
Raises:
ValueError – For unknown source or destination coordinate system or a missing
or badly formed latitude or datetime input
precision (float, optional) – Precision of output (degrees). A negative value of this argument
produces a low-precision calculation of geodetic lat/lon based only
on their spherical harmonic representation. A positive value causes
the underlying Fortran routine to iterate until feeding the output
geo lat/lon into geo2qd (APXG2Q) reproduces the input QD lat/lon to
within the specified precision.
Returns:
glat (ndarray or float) – Geodetic latitude
glon (ndarray or float) – Geodetic longitude
error (ndarray or float) – The angular difference (degrees) between the input QD coordinates
and the qlat/qlon produced by feeding the output glat and glon
into geo2qd (APXG2Q)
precision (float, optional) – Precision of output (degrees). A negative value of this argument
produces a low-precision calculation of geodetic lat/lon based only
on their spherical harmonic representation. A positive value causes
the underlying Fortran routine to iterate until feeding the output
geo lat/lon into geo2qd (APXG2Q) reproduces the input QD lat/lon to
within the specified precision.
Returns:
glat (ndarray or float) – Geodetic latitude
glon (ndarray or float) – Geodetic longitude
error (ndarray or float) – The angular difference (degrees) between the input QD coordinates
and the qlat/qlon produced by feeding the output glat and glon
into geo2qd (APXG2Q)
Computes the magnetic local time at the specified magnetic longitude
and UT.
Parameters:
mlon (array_like) – Magnetic longitude (apex and quasi-dipole longitude are always
equal)
dtime (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.
The current default is 50 * 6371, roughly 50 RE. (default=318550)
Returns:
mlt (ndarray or float) – Magnetic local time in hours [0, 24)
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).
Computes the magnetic longitude at the specified MLT and UT.
Parameters:
mlt (array_like) – Magnetic local time
dtime (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.
The current default is 50 * 6371, roughly 50 RE. (default=318550)
Returns:
mlon (ndarray or float) – Magnetic longitude [0, 360) (apex and quasi-dipole longitude are
always equal)
Notes
To compute the magnetic longitude, we find the apex longitude of the
subsolar point at the given time. Then the magnetic longitude of the
given point will be computed from the separation in magnetic local time
from this point (1 hour = 15 degrees).
Performs mapping of points along the magnetic field to the closest
or conjugate hemisphere.
Parameters:
glat (array_like) – Geodetic latitude
glon (array_like) – Geodetic longitude
height (array_like) – Source altitude in km
newheight (array_like) – Destination altitude in km
conjugate (bool, optional) – Map to newheight in the conjugate hemisphere instead of the
closest hemisphere
precision (float, optional) – Precision of output (degrees). A negative value of this argument
produces a low-precision calculation of geodetic lat/lon based only
on their spherical harmonic representation. A positive value causes
the underlying Fortran routine to iterate until feeding the output
geo lat/lon into geo2qd (APXG2Q) reproduces the input QD lat/lon to
within the specified precision.
Returns:
newglat (ndarray or float) – Geodetic latitude of mapped point
newglon (ndarray or float) – Geodetic longitude of mapped point
error (ndarray or float) – The angular difference (degrees) between the input QD coordinates
and the qlat/qlon produced by feeding the output glat and glon
into geo2qd (APXG2Q)
Notes
The mapping is done by converting glat/glon/height to modified apex
lat/lon, and converting back to geographic using newheight (if
conjugate, use negative apex latitude when converting back)
Get quasi-dipole base vectors f1 and f2 at the specified coordinates.
Parameters:
lat ((N,) array_like or float) – Latitude
lon ((N,) array_like or float) – Longitude
height ((N,) array_like or float) – Altitude in km
coords ({‘geo’, ‘apex’, ‘qd’}, optional) – Input coordinate system
precision (float, optional) – Precision of output (degrees) when converting to geo. A negative
value of this argument produces a low-precision calculation of
geodetic lat/lon based only on their spherical harmonic
representation.
A positive value causes the underlying Fortran routine to iterate
until feeding the output geo lat/lon into geo2qd (APXG2Q) reproduces
the input QD lat/lon to within the specified precision (all
coordinates being converted to geo are converted to QD first and
passed through APXG2Q).
Returns:
f1 ((2, N) or (2,) ndarray)
f2 ((2, N) or (2,) ndarray)
Notes
The vectors are described by Richmond [1995] [2] and
Emmert et al. [2010] [3]. The vector components are geodetic east and
north.
Returns base vectors in quasi-dipole and apex coordinates.
Parameters:
lat (array_like or float) – Latitude in degrees; input must be broadcastable with lon and
height.
lon (array_like or float) – Longitude in degrees; input must be broadcastable with lat and
height.
height (array_like or float) – Altitude in km; input must be broadcastable with lon and lat.
coords (str, optional) – Input coordinate system, expects one of ‘geo’, ‘apex’, or ‘qd’
(default=’geo’)
precision (float, optional) – Precision of output (degrees) when converting to geo. A negative
value of this argument produces a low-precision calculation of
geodetic lat/lon based only on their spherical harmonic
representation.
A positive value causes the underlying Fortran routine to iterate
until feeding the output geo lat/lon into geo2qd (APXG2Q) reproduces
the input QD lat/lon to within the specified precision (all
coordinates being converted to geo are converted to QD first and
passed through APXG2Q).
Returns:
f1 ((3, N) or (3,) ndarray) – Quasi-dipole base vector equivalent to e1, tangent to contours of
constant lambdaA
f2 ((3, N) or (3,) ndarray) – Quasi-dipole base vector equivalent to e2
f3 ((3, N) or (3,) ndarray) – Quasi-dipole base vector equivalent to e3, tangent to contours of
PhiA
g1 ((3, N) or (3,) ndarray) – Quasi-dipole base vector equivalent to d1
g2 ((3, N) or (3,) ndarray) – Quasi-dipole base vector equivalent to d2
g3 ((3, N) or (3,) ndarray) – Quasi-dipole base vector equivalent to d3
d1 ((3, N) or (3,) ndarray) – Apex base vector normal to contours of constant PhiA
d2 ((3, N) or (3,) ndarray) – Apex base vector that completes the right-handed system
d3 ((3, N) or (3,) ndarray) – Apex base vector normal to contours of constant lambdaA
e1 ((3, N) or (3,) ndarray) – Apex base vector tangent to contours of constant V0
e2 ((3, N) or (3,) ndarray) – Apex base vector that completes the right-handed system
e3 ((3, N) or (3,) ndarray) – Apex base vector tangent to contours of constant PhiA
Notes
The vectors are described by Richmond [1995] [4] and
Emmert et al. [2010] [5]. The vector components are geodetic east,
north, and up (only east and north for f1 and f2).
f3, g1, g2, and g3 are not part of the Fortran code
by Emmert et al. [2010] [5]. They are calculated by this
Python library according to the following equations in
Richmond [1995] [4]:
Returns the magnetic field vectors in apex coordinates.
Parameters:
lat ((N,) array_like or float) – Latitude
lon ((N,) array_like or float) – Longitude
height ((N,) array_like or float) – Altitude in km
coords ({‘geo’, ‘apex’, ‘qd’}, optional) – Input coordinate system
precision (float, optional) – Precision of output (degrees) when converting to geo. A negative
value of this argument produces a low-precision calculation of
geodetic lat/lon based only on their spherical harmonic
representation.
A positive value causes the underlying Fortran routine to iterate
until feeding the output geo lat/lon into geo2qd (APXG2Q) reproduces
the input QD lat/lon to within the specified precision (all
coordinates being converted to geo are converted to QD first and
passed through APXG2Q).
Returns:
main_mag_e3 ((1, N) or (1,) ndarray) – IGRF magnitude divided by a scaling factor, D (d_scale) to give
the main B field magnitude along the e3 base vector
e3 ((3, N) or (3,) ndarray) – Base vector tangent to the contours of constant V_0 and Phi_A
main_mag_d3 ((1, N) or (1,) ndarray) – IGRF magnitude multiplied by a scaling factor, D (d_scale) to give
the main B field magnitudee along the d3 base vector
d3 ((3, N) or (3,) ndarray) – Base vector equivalent to the scaled main field unit vector
The apex magnetic field vectors described by Richmond [1995] [4] and
Emmert et al. [2010] [5], specfically the Be3 (main_mag_e3) and Bd3
(main_mag_d3) components. The vector components are geodetic east,
north, and up.
References
Richmond, A. D. (1995) [4]
Emmert, J. T. et al. (2010) [5]