Coordinate Systems and Transformations
This package defines common coordinate systems used in heliophysics and space physics research.
Standard Coordinate Systems
Systems based on the Earth-Sun line
- GSE (Geocentric Solar Ecliptic)
- GSM (Geocentric Solar Magnetic)
Systems based on the Earth's rotation axis
- GEO (Geographic)
- GEI (Geocentric Equatorial Inertial)
- J2000
Systems based on the dipole axis of the Earth's magnetic field
- SM (Solar Magnetic)
- MAG (Geomagnetic)
Other coordinate systems
More information can be found in the the following links
Coordinate Transformations
SPEDAS.rotate
— Functionrotate(ts::AbstractMatrix, mat::AbstractMatrix)
Coordinate-aware transformation of vector/matrix by rotation matrix(s) mat(s)
. Assume ts
is a matrix of shape (n, 3).
A comprehensive description of the transformations can be found in Hapgood [1]
Coordinate transformations between geocentric systems
SPEDAS.cotrans
— Functioncotrans(A, in, out)
cotrans(A, out; in=get_coord(A))
Transform the data to the out
coordinate system from the in
coordinate system.
This function automatically choose between Julia's GeoCotrans
(if available) and Fortran's IRBEM
implementation.
References:
- IRBEM-LIB: compute magnetic coordinates and perform coordinate conversions (Documentation, IRBEM.jl)
- SPEDAS Cotrans
GeoCotrans.GeoCotrans
— ModuleCoordinate systems, transformations, and geomagnetic field models.
Coordinate systems
GEO
: Geocentric geographic cartesian coordinate system(x [𝐋], y [𝐋], z [𝐋])
.GSM
: Geocentric Solar Magnetic (GSM) coordinate system.GSE
: Geocentric Solar Ecliptic (GSE) coordinate system.GEI
: Geocentric Equatorial Inertial (GEI) coordinate system.GDZ
: Geodetic (GDZ) coordinate system(altitude [km], latitude [deg], longitude [deg])
.MAG
: Magnetic (MAG) coordinate system.SPH
: Geocentric geographic spherical (SPH) coordinate system(r [𝐋], θ [deg], φ [deg])
.J2000
: J2000 (J2000) coordinate system.
Coordinate transformations
geo2gei
,gei2geo
: Transform between GEO and GEI coordinate systems.geo2gsm
,gsm2geo
: Transform between GEO and GSM coordinate systems.gei2gsm
,gsm2gei
: Transform between GEI and GSM coordinate systems.gse2gsm
,gsm2gse
: Transform between GSE and GSM coordinate systems.
References
International Geomagnetic Reference Field (IGRF)
The International Geomagnetic Reference Field (IGRF) is a standard mathematical description of the Earth's main magnetic field. It is used widely in studies of the Earth's deep interior, crust, ionosphere, and magnetosphere.
Functions
igrf_B
: Compute the geomagnetic field (IGRF-14, dipole model)
Examples
r, θ, φ = 6500., 30., 4.
t = Date(2021, 3, 28)
Br, Bθ, Bφ = igrf_B(r, θ, φ, t)
# Input position in geodetic coordinates, output magnetic field in East-North-Up (ENU) coordinates
Be, Bn, Bu = igrf_B(GDZ(0, 60.39299, 5.32415), t)
References
Elsewhere
- SatelliteToolboxGeomagneticField.jl: Models to compute the geomagnetic field (IGRF-13, dipole model)
- ppigrf: Pure Python code to calculate IGRF model predictions.
- geopack: Python code to calculate IGRF model predictions.
using DimensionalData
using Speasy, SPEDAS
using CairoMakie
pos_gse = get_data("cda/THC_L1_STATE/thc_pos_gse", "2015-10-16", "2015-10-17") |> DimArray
pos_gsm = cotrans(pos_gse, "GSM")
pos_sm = cotrans(pos_gse, "SM")
pos_geo = cotrans(pos_gse, "GEO")
tplot((pos_gse, pos_gsm, pos_sm, pos_geo))

Specialized Coordinate Systems
The package also provides transformations for analysis-specific coordinate systems:
Field-Aligned Coordinates (FAC)
A local coordinate system defined relative to the ambient magnetic field direction, useful for studying plasma waves and particle distributions.
SPEDAS.fac_mat
— Functionfac_mat(vec::AbstractVector; xref=[1.0, 0.0, 0.0])
Generates a field-aligned coordinate (FAC) transformation matrix for a vector.
Arguments
vec
: A 3-element vector representing the magnetic field
Minimum Variance Analysis (MVA) and Boundary Normal Coordinates (LMN)
A coordinate system derived from the eigenvalues and eigenvectors of the magnetic field variance matrix, commonly used in analyzing current sheets, discontinuities, and wave propagation.
References:
- Minimum and Maximum Variance Analysis
- https://pyspedas.readthedocs.io/en/latest/coords.html#pyspedas.minvar
SPEDAS.mva_eigen
— Functionmva_eigen(B::AbstractMatrix; sort=(;), check=false) -> F::Eigen
Perform minimum variance analysis, returning Eigen
factorization object F
which contains the eigenvalues in F.values
and the eigenvectors in the columns of the matrix F.vectors
.
Set check=true
to check the reliability of the result.
The k
th eigenvector can be obtained from the slice F.vectors[:, k]
.
SPEDAS.mva
— Functionmva(V, B=V; kwargs...)
Rotate a timeseries V
into the LMN coordinates based on the reference field B
.
Arguments
V
: The timeseries data to be transformed, where each column represents a componentB
: The reference field used to determine the minimum variance directions, where each column represents a component
SPEDAS.check_mva_eigen
— Functioncheck_mva_eigen(F; r=5, verbose=false)
Check the quality of the MVA result.
If λ₁ ≥ λ₂ ≥ λ₃ are 3 eigenvalues of the constructed matrix M, then a good indicator of nice fitting LMN coordinate system should have λ₂ / λ₃ > r.
Error estimates for MVA:
SPEDAS.Δφij
— FunctionΔφij(λᵢ, λⱼ, λ₃, M)
Calculate the phase error between components i and j according to: |Δφᵢⱼ| = |Δφⱼᵢ| = √(λ₃/(M-1) * (λᵢ + λⱼ - λ₃)/(λᵢ - λⱼ)²)
Parameters:
- λᵢ: eigenvalue i
- λⱼ: eigenvalue j
- λ₃: smallest eigenvalue (λ₃)
- M: number of samples
SPEDAS.B_x3_error
— FunctionCalculate the composite statistical error estimate for ⟨B·x₃⟩: |Δ⟨B·x₃⟩| = √(λ₃/(M-1) + (Δφ₃₂⟨B⟩·x₂)² + (Δφ₃₁⟨B⟩·x₁)²)
Parameters:
- λ₁, λ₂, λ₃: eigenvalues in descending order
- M: number of samples
- B: mean magnetic field vector
- x₁, x₂, x₃: eigenvectors
See also: Comparison with PySPEDAS.
Reference
- PySPEDAS: Coordinate Systems
- geopack: Python version of geopack and Tsyganenko models
- geospacelab: A python-based library to collect, manage, and visualize geospace data (e.g. OMNI, geomagnetic indices, EISCAT, DMSP, SWARM, TEC, AMPERE, etc.).
- aacgmv2: Python library for AACGM-v2 magnetic coordinates