ELFINData

DOI

ELFINData.jl provides a high-level Julia interface to the ELFIN mission's particle and field measurements. The sections below highlight the most common entry points and link to the auto-generated API reference.

ELFINData.ELFINDataModule

Load and process data from the Electron Losses and Fields INvestigation (ELFIN) mission.

Instruments

Functions

  • epd_spectral: Load ELFIN EPD L2 data and extract directionally resolved flux spectra (omni, para, anti) and/or pitch angle spectra.

References: Website, NASA Science, Wikipedia, DOI

source

Installation

using Pkg
Pkg.add("ELFINData")

Quick Start

The examples in this section walk through the typical workflow: discover an instrument's datasets, load data for the specified time range, and request processed data products.

using ELFINData

# Inspect the Energetic Particle Detector (EPD) datasets
EPD.datasets
Dict{Tuple{String, String, String}, ELFINData.ELFINLogicalDataset{Symbol}} with 6 entries:
  ("a", "l2", "epdef") => ELA_L2_EPDEF
  ("b", "l1", "epdif") => ELB_L1_EPDIF
  ("b", "l1", "epdef") => ELB_L1_EPDEF
  ("b", "l2", "epdef") => ELB_L2_EPDEF
  ("a", "l1", "epdif") => ELA_L1_EPDIF
  ("a", "l1", "epdef") => ELA_L1_EPDEF

From the Julia REPL you can also type ?EPD to view the EPD documentation.

# Load the instrument's datasets (using the logical source name is the most portable approach)
ELA_L1_FGS("2020-10-01", "2020-10-02")
# Alternatively, specify probe and datatype explicitly. This resolves the same logical dataset.
# If no time range is given, the call returns the set of datasets (i.e. `FGM(; probe = "a", datatype = "survey") == ELA_L1_FGS`)
FGM("2020-10-01", "2020-10-02"; probe = "a", datatype = "survey")
Dataset: 
Group: ela_l1_fgs

Variables
  ela_fgs_time   (28090)
    Datatype:    CommonDataFormat.TT2000
    Dimensions:  
    Attributes:
     FILLVAL              = CommonDataFormat.TT2000[1707-09-22T12:12:10.961]
     FIELDNAM             = Time since Jan 1, 2000
     VALIDMAX             = CommonDataFormat.TT2000[2029-12-31T23:59:58.999]
     CATDESC              = ela_fgs_time, UTC, in seconds
     VALIDMIN             = CommonDataFormat.TT2000[2010-01-01T00:00:00]
     UNITS                = s
     VAR_NOTES            = Unleaped seconds
     SCALETYP             = linear
     TIME_BASE            = J2000
     FORMAT               = E12.2
     VAR_TYPE             = support_data
     LABLAXIS             = ela_fgs_time

  ela_fgs   (3 × 28090)
    Datatype:    Float32
    Dimensions:  nothing × ela_fgs_time
    Attributes:
     FILLVAL              = Float32[NaN]
     FIELDNAM             = Bfield XYZ
     VALIDMAX             = Float32[1.0f6, 1.0f6, 1.0f6]
     CATDESC              = Magnetic field B in XYZ Sensor Coordinates, Survey Mode
     DISPLAY_TYPE         = time_series
     VALIDMIN             = Float32[-1.0f6, -1.0f6, -1.0f6]
     UNITS                = nT
     VAR_NOTES            = FGM data, sensor coordinates, 10Hz
     DEPEND_0             = ela_fgs_time
     LABL_PTR_1           = StaticStrings.StaticString{6}["Bx    ", "By    ", "Bz    "]
     FORMAT               = F12.8
     VAR_TYPE             = data

  ela_fgs_labl   (3 × 1)
    Datatype:    StaticStrings.StaticString{6}
    Dimensions:  
    Attributes:
     FORMAT               = A6
     VAR_TYPE             = metadata
     CATDESC              = labels for components of Bfield
     FIELDNAM             = label Bfield
     LABLAXIS             = ela_fgs_labl

  ela_fgs_fsp_time   (945)
    Datatype:    CommonDataFormat.TT2000
    Dimensions:  
    Attributes:
     FILLVAL              = CommonDataFormat.TT2000[1707-09-22T12:12:10.961]
     FIELDNAM             = Time since Jan 1, 2000
     VALIDMAX             = CommonDataFormat.TT2000[2029-12-31T23:59:58.999]
     CATDESC              = ela_fgs_fap time, UTC, in seconds
     VALIDMIN             = CommonDataFormat.TT2000[2010-01-01T00:00:00]
     UNITS                = s
     VAR_NOTES            = Unleaped seconds
     SCALETYP             = linear
     TIME_BASE            = J2000
     FORMAT               = E12.2
     VAR_TYPE             = support_data
     LABLAXIS             = ela_fgs_fsp_time

  ela_fgs_fsp_res_dmxl   (3 × 945)
    Datatype:    Float32
    Dimensions:  nothing × ela_fgs_fsp_time
    Attributes:
     FILLVAL              = Float32[NaN]
     FIELDNAM             = FGS fsp XYZ
     VALIDMAX             = Float32[1.0f6, 1.0f6, 1.0f6]
     CATDESC              = Difference of calibrated FGM field from  IGRF fields, spin resolution
     DISPLAY_TYPE         = time_series
     VALIDMIN             = Float32[-1.0f6, -1.0f6, -1.0f6]
     UNITS                = nT
     VAR_NOTES            = FGM data, spin resoultion
     DEPEND_0             = ela_fgs_fsp_time
     LABL_PTR_1           = StaticStrings.StaticString{6}["X     ", "Y     ", "Z     "]
     COORDINATE_SYSTEM    = DMXL
     FORMAT               = F12.8
     VAR_TYPE             = data

  ela_fgs_fsp_res_dmxl_trend   (3 × 945)
    Datatype:    Float32
    Dimensions:  nothing × ela_fgs_fsp_time
    Attributes:
     FILLVAL              = Float32[NaN]
     FIELDNAM             = Detrended FGS fsp XYZ
     VALIDMAX             = Float32[1.0f6, 1.0f6, 1.0f6]
     CATDESC              = Detrended ela_fgs_fsp res_dmxl spin resolution
     DISPLAY_TYPE         = time_series
     VALIDMIN             = Float32[-1.0f6, -1.0f6, -1.0f6]
     UNITS                = nT
     VAR_NOTES            = Detrended fgs_fsp_res_dmxl data, spin resoultion
     DEPEND_0             = ela_fgs_fsp_time
     LABL_PTR_1           = StaticStrings.StaticString{6}["X     ", "Y     ", "Z     "]
     COORDINATE_SYSTEM    = DMXL
     FORMAT               = F12.8
     VAR_TYPE             = data

  ela_fgs_fsp_res_gei   (3 × 945)
    Datatype:    Float32
    Dimensions:  nothing × ela_fgs_fsp_time
    Attributes:
     FILLVAL              = Float32[NaN]
     FIELDNAM             = FGS fsp XYZ
     VALIDMAX             = Float32[1.0f6, 1.0f6, 1.0f6]
     CATDESC              = Difference of calibrated FGM field from  IGRF fields, spin resolution
     DISPLAY_TYPE         = time_series
     VALIDMIN             = Float32[-1.0f6, -1.0f6, -1.0f6]
     UNITS                = nT
     VAR_NOTES            = FGM data, spin resoultion
     DEPEND_0             = ela_fgs_fsp_time
     LABL_PTR_1           = StaticStrings.StaticString{6}["X     ", "Y     ", "Z     "]
     COORDINATE_SYSTEM    = GEI
     FORMAT               = F12.8
     VAR_TYPE             = data

  ela_fgs_fsp_labl   (3 × 1)
    Datatype:    StaticStrings.StaticString{6}
    Dimensions:  
    Attributes:
     FORMAT               = A6
     VAR_TYPE             = metadata
     CATDESC              = labels for components of fgs_fsp
     FIELDNAM             = label fgs_fsp XYZ
     LABLAXIS             = ela_fgs_fsp_labl

  ela_fgs_fsp_igrf_dmxl   (3 × 945)
    Datatype:    Float32
    Dimensions:  nothing × ela_fgs_fsp_time
    Attributes:
     FILLVAL              = Float32[NaN]
     FIELDNAM             = IGRF XYZ
     VALIDMAX             = Float32[1.0f6, 1.0f6, 1.0f6]
     CATDESC              = IGRF fields in DMXL coordinates
     DISPLAY_TYPE         = time_series
     VALIDMIN             = Float32[-1.0f6, -1.0f6, -1.0f6]
     UNITS                = nT
     VAR_NOTES            = IGRF, sampled at spin resoultion
     DEPEND_0             = ela_fgs_fsp_time
     LABL_PTR_1           = StaticStrings.StaticString{6}["X     ", "Y     ", "Z     "]
     COORDINATE_SYSTEM    = DMXL
     FORMAT               = F12.8
     VAR_TYPE             = data

  ela_fgs_fsp_igrf_gei   (3 × 945)
    Datatype:    Float32
    Dimensions:  nothing × ela_fgs_fsp_time
    Attributes:
     FILLVAL              = Float32[NaN]
     FIELDNAM             = IGRF XYZ
     VALIDMAX             = Float32[1.0f6, 1.0f6, 1.0f6]
     CATDESC              = IGRF fields in GEI coordinates
     DISPLAY_TYPE         = time_series
     VALIDMIN             = Float32[-1.0f6, -1.0f6, -1.0f6]
     UNITS                = nT
     VAR_NOTES            = IGRF, sampled at spin resoultion
     DEPEND_0             = ela_fgs_fsp_time
     LABL_PTR_1           = StaticStrings.StaticString{6}["X     ", "Y     ", "Z     "]
     COORDINATE_SYSTEM    = GEI
     FORMAT               = F12.8
     VAR_TYPE             = data

  ela_fgs_igrf_labl   (3 × 1)
    Datatype:    StaticStrings.StaticString{6}
    Dimensions:  
    Attributes:
     FORMAT               = A6
     VAR_TYPE             = metadata
     CATDESC              = labels for components of igrf
     FIELDNAM             = label fsp XYZ
     LABLAXIS             = ela_fgs_igrf_labl

  ela_fgs_att_flag   (3 × 0)
    Datatype:    Int8
    Dimensions:  nothing × ela_fgs_fsp_time
    Attributes:
     FILLVAL              = Int8[0]
     FIELDNAM             = Attitude flag
     VALIDMAX             = Int8[1]
     CATDESC              = Flag for attitude model1ed vs interpolated - interpolated: 0, modelled: 1
     DISPLAY_TYPE         = time_series
     VALIDMIN             = Int8[0]
     UNITS                =  
     VAR_NOTES            = Attitude flag,  modelled or interpolated
     DEPEND_0             = ela_fgs_fsp_time
     FORMAT               = i1
     VAR_TYPE             = data

Global attributes
  Project              = ["Electron Losses and Fields INvestigation"]
  Source_name          = ["ela>Electron Losses and Fields INvestigation-A"]
  Discipline           = ["Space Physics>Magnetospheric Science", "Space Physics>Interplanetary Studies"]
  Data_type            = ["FGM"]
  Descriptor           = ["L1>L1 DATA"]
  File_naming_convention = ["source_descriptor_datatype"]
  Data_version         = ["1"]
  PI_name              = ["V. Angelopoulos"]
  PI_affiliation       = ["UCLA, IGPP/EPSS"]
  TEXT                 = ["The Electron Losses and Fields Investigation (ELFIN) mission is a space weather mission using three scientific instruments in a 3U+ CubeSat. The instruments measure wave and particle data."]
  Instrument_type      = ["Magnetic Fields (space)"]
  Mission_group        = ["ELFIN"]
  Logical_source       = ["ela_l1_fgs"]
  Logical_file_id      = ["ela_l1_fgs_00000000_v01"]
  Logical_source_description = ["Spacecraft fluxgate magnetometer, Survey mode, raw sensor data"]
  Time_resolution      = ["10hz"]
  Rules_of_use         = ["Open Data for Scientific Use"]
  Generated_by         = ["ELFIN SOC"]
  Generation_date      = ["2018-10-20"]
  Acknowledgement      = ["NSF, NASA "]
  MODS                 = ["Rev- 2018-10-20"]
  ADID_ref             = ["xxx"]
  LINK_TEXT            = ["http://elfin.igpp.ucla.edu"]
  LINK_TITLE           = ["Electron Losses and Fields INvestigation-A FGM sensor data"]
  HTTP_LINK            = ["http://elfin.igpp.ucla.edu"]

The following function derives directionally resolved flux spectra (omni, para, anti) and/or pitch-angle spectra from EPD level 2 data with a single call.

epd_spectral("2020-10-01", "2020-10-02"; probe = "a")
16×1731 DimStack
├──────────────────┴───────────────────────────────────────────────────── dims ┐
  ↓ Energy Sampled{Float32} [63.24554f0, …, 6500.0f0] ForwardOrdered Irregular Points,
  → Ti Sampled{Dates.DateTime} [Dates.DateTime("2020-10-01T00:26:53.194"), …, Dates.DateTime("2020-10-01T20:45:49.320")] ForwardOrdered Irregular Points
├────────────────────────────────────────────────────────────────────── layers ┤
  :omni eltype: Float32 dims: Energy, Ti size: 16×1731
  :para eltype: Float32 dims: Energy, Ti size: 16×1731
  :anti eltype: Float32 dims: Energy, Ti size: 16×1731
  :perp eltype: Float32 dims: Energy, Ti size: 16×1731
  :prec eltype: Float32 dims: Energy, Ti size: 16×1731
├──────────────────────────────────────────────────────────────────── metadata ┤
  Dict{Union{AbstractString, Symbol}, Any} with 13 entries:
  "FILLVAL"      => Float32[NaN]
  "DEPEND_0"     => "ela_pef_hs_time"
  "FIELDNAM"     => "nflux"
  "VALIDMAX"     => Float32[1.0f6]
  "DEPEND_1"     => "ela_pef_hs_epa_spec"
  "FORMAT"       => "F12.8"
  "VAR_TYPE"     => "data"
  "CATDESC"      => "half spin res, EnergyPitchAngleTime spectra nflux"
  "LABLAXIS"     => "nflux"
  "VALIDMIN"     => Float32[0.0]
  "DISPLAY_TYPE" => "spectral"
  "DEPEND_2"     => "ela_pef_energies_mean"
  "UNITS"        => "#/(cm^2*s*str*MeV)"
└──────────────────────────────────────────────────────────────────────────────┘

Commonly used variables have concise convenience wrappers. Because these variables are uniquely named, no additional metadata (such as probe, datatype, level, or dataset) is required to access them.

ELA_POS_GEI("2020-10-01", "2020-10-02")
# Alternative ways to access the same variable:
# `STATE("2020-10-01", "2020-10-02"; probe = "a")["ela_pos_gei"]` or
# `ELA_L1_STATE("2020-10-01", "2020-10-02")["ela_pos_gei"]`

ELB_FGS("2020-10-01", "2020-10-02")
# Alternative ways to access the same variable:
# `ELB_L1_FGS("2020-10-01", "2020-10-02")["elb_fgs"]` or
# `FGM("2020-10-01", "2020-10-02"; probe = "b", datatype = "survey")["elb_fgs"]`

ELA_PEF_HS_EPAT_NFLUX("2020-10-01", "2020-10-02")
# Alternative ways to access the same variable:
# `ELA_L2_EPDEF("2020-10-01", "2020-10-02")["ela_pef_hs_Epat_nflux"]` or
# `EPD("2020-10-01", "2020-10-02"; probe = "a", level = "l2", datatype = "epdef")["ela_pef_hs_Epat_nflux"]`
ela_pef_hs_Epat_nflux (10 × 16 × 1731)
  Datatype:    Float32
  Dimensions:  ela_pef_hs_epa_spec × ela_pef_energies_mean × ela_pef_hs_time
  Attributes:
   FILLVAL              = Float32[NaN]
   FIELDNAM             = nflux
   VALIDMAX             = Float32[1.0f6]
   DEPEND_1             = ela_pef_hs_epa_spec
   CATDESC              = half spin res, EnergyPitchAngleTime spectra nflux
   VALIDMIN             = Float32[0.0]
   DISPLAY_TYPE         = spectral
   DEPEND_2             = ela_pef_energies_mean
   UNITS                = #/(cm^2*s*str*MeV)
   DEPEND_0             = ela_pef_hs_time
   FORMAT               = F12.8
   VAR_TYPE             = data
   LABLAXIS             = nflux

API

Instruments

Datasets

Variables

Functions and Types

ELFINData.epd_spectralMethod
epd_spectral(trange; probe = "a", type = "nflux", datatype = "epdef", fullspin = false)

Load ELFIN EPD L2 data and process it to extract directionally resolved flux spectra (omni, para, anti) and/or pitch angle spectra.

source

Validation and Benchmarking with PySPEDAS

We compare ELFINData results against PySPEDAS to ensure numerical parity and to track performance. The following benchmarks illustrate the typical workflow.

State variable ELA_POS_GEI: Julia is about 200 times faster than Python for this retrieval.

using PySPEDAS
using ELFINData
using Chairmarks

trange = ("2021-08-08", "2021-08-10")
ela_pos_gei = ELA_POS_GEI(trange)
py_ela_pos_gei = Array(PySPEDAS.elfin.state(trange).ela_pos_gei)
@assert ela_pos_gei == py_ela_pos_gei'
@b Array(ELA_POS_GEI(trange)), PySPEDAS.elfin.state(trange), pyspedas.projects.elfin.state(trange)
(1.116 ms (249 allocs: 3.966 MiB, without a warmup), 3.133 s (841 allocs: 31.656 KiB, without a warmup), 3.169 s (52 allocs: 2.391 KiB, 0.14% compile time, without a warmup))

Processing EPD L2 spectra: Julia is about 100 times faster than Python for spectral derivations across the same interval.

trange = ("2020-10-01", "2020-10-02")
nflux_para = epd_spectral(trange).para
py_nflux_para = PySPEDAS.elfin.epd(trange; level = "l2").ela_pef_hs_nflux_para
@assert nflux_para ≈ Array(py_nflux_para)'
b1 = @b epd_spectral($trange)
b2 = @b PySPEDAS.elfin.epd($trange; level = "l2")
@info "Julia" b1
@info "PySPEDAS" b2
03-Nov-25 17:15:30: ELFIN EPD: START LOADING.
03-Nov-25 17:15:30: Downloading remote index: https://data.elfin.ucla.edu/ela/l2/epd/fast/electron/2020/
03-Nov-25 17:15:32: Downloading https://data.elfin.ucla.edu/ela/l2/epd/fast/electron/2020/ela_l2_epdef_20201001_v01.cdf to elfin_data/ela/l2/epd/fast/electron/2020/ela_l2_epdef_20201001_v01.cdf
03-Nov-25 17:15:32: Download of elfin_data/ela/l2/epd/fast/electron/2020/ela_l2_epdef_20201001_v01.cdf complete, 9.063 MB in 0.3 sec (29.375 MB/sec) (transfer_normal)
03-Nov-25 17:15:32: Unable to get ydata for variable ela_pef_Tspin
03-Nov-25 17:15:32: Variable ela_pef_hs_epa_spec DEPEND_1 attribute ela_pef_energies_mean has length 16, but corresponding data dimension has length 10. Removing attribute.
03-Nov-25 17:15:32: Variable ela_pef_fs_epa_spec DEPEND_1 attribute ela_pef_energies_mean has length 16, but corresponding data dimension has length 18. Removing attribute.
03-Nov-25 17:15:32: ELFIN EPD: LOADING END.
03-Nov-25 17:15:32: ELFIN EPD L2: START PROCESSING.
03-Nov-25 17:15:32: ELFIN EPD L2: START ENERGY SPECTOGRAM.
03-Nov-25 17:15:32: ELFIN EPD L2: START PITCH ANGLE SPECTOGRAM.
03-Nov-25 17:15:32: Energy channel [(0, 2), (3, 5), (6, 8), (9, 15)] are used for epd l2 pitch angle spectra.
03-Nov-25 17:15:34: ELFIN EPD: START LOADING.
03-Nov-25 17:15:34: Downloading remote index: https://data.elfin.ucla.edu/ela/l2/epd/fast/electron/2020/
03-Nov-25 17:15:35: Downloading https://data.elfin.ucla.edu/ela/l2/epd/fast/electron/2020/ela_l2_epdef_20201001_v01.cdf to elfin_data/ela/l2/epd/fast/electron/2020/ela_l2_epdef_20201001_v01.cdf
03-Nov-25 17:15:35: Download of elfin_data/ela/l2/epd/fast/electron/2020/ela_l2_epdef_20201001_v01.cdf complete, 9.063 MB in 0.3 sec (29.400 MB/sec) (transfer_normal)
03-Nov-25 17:15:36: Unable to get ydata for variable ela_pef_Tspin
03-Nov-25 17:15:36: Variable ela_pef_hs_epa_spec DEPEND_1 attribute ela_pef_energies_mean has length 16, but corresponding data dimension has length 10. Removing attribute.
03-Nov-25 17:15:36: Variable ela_pef_fs_epa_spec DEPEND_1 attribute ela_pef_energies_mean has length 16, but corresponding data dimension has length 18. Removing attribute.
03-Nov-25 17:15:36: ELFIN EPD: LOADING END.
03-Nov-25 17:15:36: ELFIN EPD L2: START PROCESSING.
03-Nov-25 17:15:36: ELFIN EPD L2: START ENERGY SPECTOGRAM.
03-Nov-25 17:15:36: ELFIN EPD L2: START PITCH ANGLE SPECTOGRAM.
03-Nov-25 17:15:36: Energy channel [(0, 2), (3, 5), (6, 8), (9, 15)] are used for epd l2 pitch angle spectra.
┌ Info: Julia
  b1 = 2.147 ms (941 allocs: 1.746 MiB)
┌ Info: PySPEDAS
  b2 = 2.043 s (1089 allocs: 44.984 KiB, without a warmup)
Array layout

Julia arrays follow the column-major convention used by most CDF files (time is the last dimension), whereas NumPy and PySPEDAS use row-major arrays (time is the first dimension). Transpose the PySPEDAS output before comparing so the dimensions align.