ELFINData
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.ELFINData — Module
Load and process data from the Electron Losses and Fields INvestigation (ELFIN) mission.
Instruments
- Fluxgate Magnetometer (FGM)
- Energetic Particle Detector (EPD)
- Magneto Resistive Magnetometer-a (MRMa)
- Magneto Resistive Magnetometer-i (MRMi)
- State data (state)
- Engineering data (ENG)
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
Installation
using Pkg
Pkg.add("ELFINData")An agent skill is included for using ELFINData with natural language. To install it using skills, run:
npx skills add JuliaSpacePhysics/ELFINData.jlQuick 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.datasetsDict{Tuple{String, String, String}, ELFINData.ELFINLogicalDataset{Symbol, MD, ELFINData.Pattern} where MD} 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_EPDEFFrom 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: ["elfin_data/ela/l1/fgm/survey/2020/ela_l1_fgs_20201001_v01.cdf"]
Group: ela_l1_fgs
Data variables
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 = CommonDataFormat.StaticString{6, UInt8}["Bx ", "By ", "Bz "]
FORMAT = F12.8
VAR_TYPE = data
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 = CommonDataFormat.StaticString{6, UInt8}["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 = CommonDataFormat.StaticString{6, UInt8}["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 = CommonDataFormat.StaticString{6, UInt8}["X ", "Y ", "Z "]
COORDINATE_SYSTEM = GEI
FORMAT = F12.8
VAR_TYPE = data
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 = CommonDataFormat.StaticString{6, UInt8}["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 = CommonDataFormat.StaticString{6, UInt8}["X ", "Y ", "Z "]
COORDINATE_SYSTEM = GEI
FORMAT = F12.8
VAR_TYPE = data
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
Support variables: ela_fgs_time, ela_fgs_fsp_time
Metadata variables: ela_fgs_labl, ela_fgs_fsp_labl, ela_fgs_igrf_labl
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{Any, Any} with 17 entries:
"SCALETYP" => log10
"FILLVAL" => Float32[NaN]
"DEPEND_0" => "ela_pef_hs_time"
"FIELDNAM" => "nflux"
"VALIDMAX" => Float32[1.0f6]
:colorrange => (10.0, 1.0e7)
"DEPEND_1" => "ela_pef_hs_epa_spec"
:yscale => log10
"FORMAT" => "F12.8"
:ylabel => "Energy (keV)"
"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
Quick Plots
using ELFINData
using Dates
using ELFINData.DimensionalData
using SpacePhysicsMakie, WGLMakie
# https://data.elfin.ucla.edu/ela/overplots/2022/09/05/ela_l2_overview_20220905_10_ndes.gif
t0 = DateTime("2022-09-05T10:00:00")
t1 = DateTime("2022-09-05T10:30:00")
# Example: plot EPD flux spectra
spectra = epd_spectral(t0, t1; probe = "a")[Ti(t0 .. t1)]
tplot(degap.([spectra.omni, spectra.anti, spectra.perp, spectra.para]); colormap=:turbo)API
Instruments
ELFINData.EPD — Constant
Energetic Particle Detector (EPD)
The EPD consists of two heads: an electron head that measures 50 keV to 5 MeV electrons and an ion head that measures 50 to 5000 keV ions.
Datasets:
- ELFIN A:
ELA_L1_EPDEF,ELA_L1_EPDIF,ELA_L2_EPDEF - ELFIN B:
ELB_L1_EPDEF,ELB_L1_EPDIF,ELB_L2_EPDEF
ELFINData.FGM — Constant
Fluxgate Magnetometer (FGM)
Datasets: ELA_L1_FGS, ELB_L1_FGS
ELFINData.STATE — Constant
State data (STATE)
Datasets: ELA_L1_STATE, ELB_L1_STATE
Datasets
ELFINData.ELA_L1_EPDEF — Constant
ELFIN A L1 EPD Electron Counts (Version 1)
ELFINData.ELA_L1_EPDIF — Constant
ELFIN A L1 EPD Ion Counts (Version 1)
ELFINData.ELA_L1_FGS — Constant
ELFIN A L1 FGM Survey (Version 1)
Main data variables: ELA_FGS
ELFINData.ELA_L1_MRMA — Constant
ELFIN A L1 MRMA (Version 1)
Main data variables: ELA_MRMA
ELFINData.ELA_L1_MRMI — Constant
ELFIN A L1 MRMI (Version 1)
Main data variables: ELA_MRMI
ELFINData.ELA_L1_STATE — Constant
ELFIN A L1 State (Version 2)
Main data variables: ELA_POS_GEI
ELFINData.ELA_L2_EPDEF — Constant
ELFIN A L2 EPD Electron Counts (Version 1)
Main data variables: ELA_PEF_HS_EPAT_NFLUX, ELA_PEF_FS_EPAT_NFLUX, ELA_PEF_HS_EPAT_EFLUX, ELA_PEF_FS_EPAT_EFLUX
ELFINData.ELB_L1_EPDEF — Constant
ELFIN B L1 EPD Electron Counts (Version 1)
ELFINData.ELB_L1_EPDIF — Constant
ELFIN B L1 EPD Ion Counts (Version 1)
ELFINData.ELB_L1_FGS — Constant
ELFIN B L1 FGM Survey (Version 1)
Main data variables: ELB_FGS
ELFINData.ELB_L1_MRMA — Constant
ELFIN B L1 MRMA (Version 1)
Main data variables: ELB_MRMA
ELFINData.ELB_L1_MRMI — Constant
ELFIN B L1 MRMI (Version 1)
Main data variables: ELB_MRMI
ELFINData.ELB_L1_STATE — Constant
ELFIN B L1 State (Version 2)
Main data variables: ELB_POS_GEI
ELFINData.ELB_L2_EPDEF — Constant
ELFIN B L2 EPD Electron Counts (Version 1)
Main data variables: ELB_PEF_HS_EPAT_NFLUX, ELB_PEF_FS_EPAT_NFLUX, ELB_PEF_HS_EPAT_EFLUX, ELB_PEF_FS_EPAT_EFLUX
Variables
ELFINData.ELA_FGS — Constant
ELFIN A FGM Magnetic field B in XYZ Sensor Coordinates, Survey Mode
ELFINData.ELA_MRMA — Constant
Magnetic field XYZ Sensor data collected by ACB, Sensor Coordinates, ADC units
ELFINData.ELA_MRMI — Constant
Magnetic field XYZ Sensor data collected by IDPU, Sensor Coordinates, ADC units
ELFINData.ELA_PEF_FS_EPAT_EFLUX — Constant
ELFIN A EnergyPitchAngleTime spectra eflux, full spin resolution
ELFINData.ELA_PEF_FS_EPAT_NFLUX — Constant
ELFIN A EnergyPitchAngleTime spectra nflux, full spin resolution
ELFINData.ELA_PEF_HS_EPAT_EFLUX — Constant
ELFIN A EnergyPitchAngleTime spectra eflux, half spin resolution
ELFINData.ELA_PEF_HS_EPAT_NFLUX — Constant
ELFIN A EnergyPitchAngleTime spectra nflux, half spin resolution
ELFINData.ELA_POS_GEI — Constant
ELFIN A State Position XYZ in GEI coordinates
ELFINData.ELB_FGS — Constant
ELFIN B FGM Magnetic field B in XYZ Sensor Coordinates, Survey Mode
ELFINData.ELB_MRMA — Constant
Magnetic field XYZ Sensor data collected by ACB, Sensor Coordinates, ADC units
ELFINData.ELB_MRMI — Constant
Magnetic field XYZ Sensor data collected by IDPU, Sensor Coordinates, ADC units
ELFINData.ELB_PEF_FS_EPAT_EFLUX — Constant
ELFIN B EnergyPitchAngleTime spectra eflux, full spin resolution
ELFINData.ELB_PEF_FS_EPAT_NFLUX — Constant
ELFIN B EnergyPitchAngleTime spectra nflux, full spin resolution
ELFINData.ELB_PEF_HS_EPAT_EFLUX — Constant
ELFIN B EnergyPitchAngleTime spectra eflux, half spin resolution
ELFINData.ELB_PEF_HS_EPAT_NFLUX — Constant
ELFIN B EnergyPitchAngleTime spectra nflux, half spin resolution
ELFINData.ELB_POS_GEI — Constant
ELFIN B State Position XYZ in GEI coordinates
Functions and Types
ELFINData.epd_spectral — Method
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.
Returns a NamedTuple-like container with omni, para, anti, and prec flux spectra.
epd_spectral("2020-10-01", "2020-10-02")