apexpy
¶
Submodules¶
Package Contents¶
Classes¶
Calculates coordinate conversions, field-line mapping, and base vectors. |
- class apexpy.Apex(date=None, refh=0, datafile=None, fortranlib=None)[source]¶
Bases:
object
Calculates coordinate conversions, field-line mapping, and base vectors.
- Parameters
date (NoneType, float,
dt.date
, ordt.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.
The geodetic reference ellipsoid is WGS84.
References
- 1
Thébault, E. et al. (2015), International Geomagnetic Reference Field: the 12th generation, Earth, Planets and Space, 67(1), 79, doi:10.1186/s40623-015-0228-9.
- __repr__()¶
Produce an evaluatable representation of the Apex class.
- __str__()¶
Produce a user-friendly representation of the Apex class.
- __eq__(comp_obj)¶
Performs equivalency evaluation.
- Parameters
comp_obj – Object of any time to be compared to the class object
- Returns
bool or NotImplemented – True if self and comp_obj are identical, False if they are not, and NotImplemented if the classes are not the same
- _apex2qd_nonvectorized(alat, alon, height)¶
Convert from apex to quasi-dipole (not-vectorised)
- Parameters
alat ((float)) – Apex latitude in degrees
alon ((float)) – Apex longitude in degrees
height ((float)) – Height in km
- Returns
qlat ((float)) – Quasi-dipole latitude in degrees
qlon ((float)) – Quasi-diplole longitude in degrees
- _qd2apex_nonvectorized(qlat, qlon, height)¶
Converts quasi-dipole to modified apex coordinates.
- Parameters
qlat (float) – Quasi-dipole latitude
qlon (float) – Quasi-dipole longitude
height (float) – Altitude in km
- Returns
alat (float) – Modified apex latitude
alon (float) – Modified apex longitude
- Raises
ApexHeightError – if apex height < reference height
- _map_EV_to_height(alat, alon, height, newheight, data, ev_flag)¶
Maps electric field related values to a desired height
- Parameters
alat (array-like) – Apex latitude in degrees.
alon (array-like) – Apex longitude in degrees.
height (array-like) – Current altitude in km.
new_height (array-like) – Desired altitude to which EV values will be mapped in km.
data (array-like) – 3D value(s) for the electric field or electric drift
ev_flag (str) – Specify if value is an electric field (‘E’) or electric drift (‘V’)
- Returns
data_mapped (array-like) – Data mapped along the magnetic field from the old height to the new height.
- Raises
ValueError – If an unknown ev_flag or badly shaped data input is supplied.
- _get_babs_nonvectorized(glat, glon, height)¶
Get the absolute value of the B-field in Tesla
- Parameters
glat (float) – Geodetic latitude in degrees
glon (float) – Geodetic longitude in degrees
height (float) – Altitude in km
- Returns
babs (float) – Absolute value of the magnetic field in Tesla
- convert(lat, lon, source, dest, height=0, datetime=None, precision=1e-10, ssheight=50 * 6371)¶
Converts between geodetic, modified apex, quasi-dipole and MLT.
- Parameters
lat (array_like) – Latitude in degrees
lon (array_like) – Longitude in degrees or MLT in hours
source (str) – Input coordinate system, accepts ‘geo’, ‘apex’, ‘qd’, or ‘mlt’
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
- geo2apex(glat, glon, height)¶
Converts geodetic to modified apex coordinates.
- Parameters
glat (array_like) – Geodetic latitude
glon (array_like) – Geodetic longitude
height (array_like) – Altitude in km
- Returns
alat (ndarray or float) – Modified apex latitude
alon (ndarray or float) – Modified apex longitude
- apex2geo(alat, alon, height, precision=1e-10)¶
Converts modified apex to geodetic coordinates.
- Parameters
alat (array_like) – Modified apex latitude
alon (array_like) – Modified apex longitude
height (array_like) – Altitude in km
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)
- geo2qd(glat, glon, height)¶
Converts geodetic to quasi-dipole coordinates.
- Parameters
glat (array_like) – Geodetic latitude
glon (array_like) – Geodetic longitude
height (array_like) – Altitude in km
- Returns
qlat (ndarray or float) – Quasi-dipole latitude
qlon (ndarray or float) – Quasi-dipole longitude
- qd2geo(qlat, qlon, height, precision=1e-10)¶
Converts quasi-dipole to geodetic coordinates.
- Parameters
qlat (array_like) – Quasi-dipole latitude
qlon (array_like) – Quasi-dipole longitude
height (array_like) – Altitude in km
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)
- apex2qd(alat, alon, height)¶
Converts modified apex to quasi-dipole coordinates.
- Parameters
alat (array_like) – Modified apex latitude
alon (array_like) – Modified apex longitude
height (array_like) – Altitude in km
- Returns
qlat (ndarray or float) – Quasi-dipole latitude
qlon (ndarray or float) – Quasi-dipole longitude
- Raises
ApexHeightError – if height > apex height
- qd2apex(qlat, qlon, height)¶
Converts quasi-dipole to modified apex coordinates.
- Parameters
qlat (array_like) – Quasi-dipole latitude
qlon (array_like) – Quasi-dipole longitude
height (array_like) – Altitude in km
- Returns
alat (ndarray or float) – Modified apex latitude
alon (ndarray or float) – Modified apex longitude
- Raises
ApexHeightError – if apex height < reference height
- mlon2mlt(mlon, dtime, ssheight=318550)¶
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 timessheight (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).
- mlt2mlon(mlt, dtime, ssheight=318550)¶
Computes the magnetic longitude at the specified MLT and UT.
- Parameters
mlt (array_like) – Magnetic local time
dtime (
datetime.datetime
) – Date and timessheight (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).
- map_to_height(glat, glon, height, newheight, conjugate=False, precision=1e-10)¶
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)
- map_E_to_height(alat, alon, height, newheight, edata)¶
Performs mapping of electric field along the magnetic field.
It is assumed that the electric field is perpendicular to B.
- Parameters
alat ((N,) array_like or float) – Modified apex latitude
alon ((N,) array_like or float) – Modified apex longitude
height ((N,) array_like or float) – Source altitude in km
newheight ((N,) array_like or float) – Destination altitude in km
edata ((3,) or (3, N) array_like) – Electric field (at alat, alon, height) in geodetic east, north, and up components
- Returns
out ((3, N) or (3,) ndarray) – The electric field at newheight (geodetic east, north, and up components)
- map_V_to_height(alat, alon, height, newheight, vdata)¶
Performs mapping of electric drift velocity along the magnetic field.
It is assumed that the electric field is perpendicular to B.
- Parameters
alat ((N,) array_like or float) – Modified apex latitude
alon ((N,) array_like or float) – Modified apex longitude
height ((N,) array_like or float) – Source altitude in km
newheight ((N,) array_like or float) – Destination altitude in km
vdata ((3,) or (3, N) array_like) – Electric drift velocity (at alat, alon, height) in geodetic east, north, and up components
- Returns
out ((3, N) or (3,) ndarray) – The electric drift velocity at newheight (geodetic east, north, and up components)
- basevectors_qd(lat, lon, height, coords='geo', precision=1e-10)¶
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.
References
- 2
Richmond, A. D. (1995), Ionospheric Electrodynamics Using Magnetic Apex Coordinates, Journal of geomagnetism and geoelectricity, 47(2), 191–212, doi:10.5636/jgg.47.191.
- 3
Emmert, J. T., A. D. Richmond, and D. P. Drob (2010), A computationally compact representation of Magnetic-Apex and Quasi-Dipole coordinates with smooth base vectors, J. Geophys. Res., 115(A8), A08322, doi:10.1029/2010JA015326.
- basevectors_apex(lat, lon, height, coords='geo', precision=1e-10)¶
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:
g1: Eqn. 6.3
g2: Eqn. 6.4
g3: Eqn. 6.5
f3: Eqn. 6.8
References
- 4(1,2,3,4,5)
Richmond, A. D. (1995), Ionospheric Electrodynamics Using Magnetic Apex Coordinates, Journal of geomagnetism and geoelectricity, 47(2), 191–212, doi:10.5636/jgg.47.191.
- 5(1,2,3,4)
Emmert, J. T., A. D. Richmond, and D. P. Drob (2010), A computationally compact representation of Magnetic-Apex and Quasi-Dipole coordinates with smooth base vectors, J. Geophys. Res., 115(A8), A08322, doi:10.1029/2010JA015326.
- get_apex(lat, height=None)¶
Calculate apex height
- Parameters
lat (float) – Apex latitude in degrees
height (float or NoneType) – Height above the surface of the Earth in km or NoneType to use reference height (default=None)
- Returns
apex_height (float) – Height of the field line apex in km
- get_height(lat, apex_height)¶
Calculate the height given an apex latitude and apex height.
- Parameters
lat (float) – Apex latitude in degrees
apex_height (float) – Maximum height of the apex field line above the surface of the Earth in km
- Returns
height (float) – Height above the surface of the Earth in km
- set_epoch(year)¶
Updates the epoch for all subsequent conversions.
- Parameters
year (float) – Decimal year
- set_refh(refh)¶
Updates the apex reference height for all subsequent conversions.
- Parameters
refh (float) – Apex reference height in km
Notes
The reference height is the height to which field lines will be mapped, and is only relevant for conversions involving apex (not quasi-dipole).
- get_babs(glat, glon, height)¶
Returns the magnitude of the IGRF magnetic field in tesla.
- Parameters
glat (array_like) – Geodetic latitude in degrees
glon (array_like) – Geodetic longitude in degrees
height (array_like) – Altitude in km
- Returns
babs (ndarray or float) – Magnitude of the IGRF magnetic field in Tesla
- bvectors_apex(lat, lon, height, coords='geo', precision=1e-10)¶
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
Notes
See Richmond, A. D. (1995) 4 equations 3.8-3.14
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