MMS Science event

Electron-scale measurements of magnetic reconnection in space

Obtain and analyze plasma data, including spectra, for the first MMS Science event (Burch et al., Science 2016, link here) showing magnetopause reconnection at the electron diffusion region on 2016/10/16 13:07:02.2 UT. Start by using Hwk06_mpause_RX.pro, provided. The objective of this exercise is to introduce plasma distributions from MMS (electrons and ions) and create burst spectra. In the process also create plasma moments and the plasma current and plot these along E and B. You are requested to create Figure 2 of Burch et al. plus the 3 bottom (electron) spectrograms from Figure 3 of the same paper (panels 3G, 3H and 3I). Your figure should look like the one shown on the next page, but using burst mode data in order to make it look like in Burch et al. You are requested to plot this figure at 2 time scales: the overview (2 min) timescale as in Fig. 2 of Burch et al. [‘13:05:30’,’13:07:30’] and the zoom-in (3 sec) timescale as in Fig. 3 of Burch et al. [‘13:07:00.5’,‘13:07:03.5’]. Notice that the clean and fast m’pause crossing was at 13:05:40UT, and this is used to determine N. In Fig. 2K of Burch et al. (the right-hand side of Fig. 2, the pictorial view of the MMS trajectory for the 2min interval) this initial m’pause crossing was near the start of the trajectory. The trajectory crosses the X-point at the 3 seconds of the zoom-in interval.

using Speasy
using DimensionalData
using SPEDAS
using SPEDAS.MMS
using SPEDAS: dimarrayify
using CairoMakie
using Unitful

update_theme!(colormap=:rainbow1)
SPEDAS.DEFAULTS.add_title = true;
    CondaPkg Found dependencies: /home/runner/.julia/packages/DimensionalData/5jhQ2/CondaPkg.toml
    CondaPkg Found dependencies: /home/runner/.julia/packages/Speasy/xrPd4/CondaPkg.toml
    CondaPkg Found dependencies: /home/runner/.julia/packages/PythonCall/L4cjh/CondaPkg.toml
    CondaPkg Initialising pixi
             /home/runner/.julia/artifacts/cefba4912c2b400756d043a2563ef77a0088866b/bin/pixi
             init
             --format pixi
             /home/runner/work/juliaspacephysics.github.io/juliaspacephysics.github.io/.CondaPkg
✔ Created /home/runner/work/juliaspacephysics.github.io/juliaspacephysics.github.io/.CondaPkg/pixi.toml
    CondaPkg Wrote /home/runner/work/juliaspacephysics.github.io/juliaspacephysics.github.io/.CondaPkg/pixi.toml
             [dependencies]
             openssl = ">=3, <3.6, >=3, <3.6"
             uv = ">=0.4"
             libstdcxx-ng = ">=3.4,<15.0"
             xarray = "*"
             sqlite = "!=3.49.1"
             numpy = "*"
                              [dependencies.python]
                 channel = "conda-forge"
                 build = "*cpython*"
                 version = ">=3.8,<4"
                          [project]
             name = ".CondaPkg"
             platforms = ["linux-64"]
             channels = ["conda-forge"]
             channel-priority = "strict"
             description = "automatically generated by CondaPkg.jl"
                          [pypi-dependencies.speasy]
             git = "https://github.com/SciQLop/speasy"
    CondaPkg Installing packages
             /home/runner/.julia/artifacts/cefba4912c2b400756d043a2563ef77a0088866b/bin/pixi
             install
             --manifest-path /home/runner/work/juliaspacephysics.github.io/juliaspacephysics.github.io/.CondaPkg/pixi.toml
✔ The default environment has been installed.
Non compliant ISTP file: No data variable found, this is suspicious
Non compliant ISTP file: No data variable found, this is suspicious
Non compliant ISTP file: No data variable found, this is suspicious
Non compliant ISTP file: No data variable found, this is suspicious
Non compliant ISTP file: No data variable found, this is suspicious
Non compliant ISTP file: Epoch was marked as data variable but it has 0 support variable

Overview plot

Minimum Variance Analysis

trange = ("2015-10-16T13:05:35", "2015-10-16T13:07:25")
tr_mpause = ("2015-10-16T13:05:40", "2015-10-16T13:06:09")

probe = 2

tvars = (
    "cda/MMS2_FGM_SRVY_L2/mms2_fgm_b_gse_srvy_l2_clean",
    "cda/MMS2_EDP_FAST_L2_DCE/mms2_edp_dce_gse_fast_l2",
)
B_gse, edp_dce_gse = get_data(tvars, trange) .|> DimArray
rotMat = mva_eigen(tview(B_gse[:, 1:3], tr_mpause))
B_LMN = rotate(B_gse[:, 1:3], rotMat) |> set_coord("LMN") |> set_coord("Boundary-normal coordinates"; old_coords=["Geocentric Solar Ecliptic"])
E_LMN = rotate(edp_dce_gse, rotMat) |> set_coord("LMN")

rotMat
Can't get MMS2_FGM_SRVY_L2/mms2_fgm_b_gse_srvy_l2_clean without web service, switching to web service
LinearAlgebra.Eigen{Float32, Quantity{Float32, 𝐌^2 𝐈^-2 𝐓^-4, Unitful.FreeUnits{(nT^2,), 𝐌^2 𝐈^-2 𝐓^-4, nothing}}, StaticArraysCore.SMatrix{3, 3, Float32, 9}, StaticArraysCore.SVector{3, Quantity{Float32, 𝐌^2 𝐈^-2 𝐓^-4, Unitful.FreeUnits{(nT^2,), 𝐌^2 𝐈^-2 𝐓^-4, nothing}}}}
values:
3-element StaticArraysCore.SVector{3, Quantity{Float32, 𝐌^2 𝐈^-2 𝐓^-4, Unitful.FreeUnits{(nT^2,), 𝐌^2 𝐈^-2 𝐓^-4, nothing}}} with indices SOneTo(3):
 501.50867f0 nT^2
  32.284077f0 nT^2
  13.253269f0 nT^2
vectors:
3×3 StaticArraysCore.SMatrix{3, 3, Float32, 9} with indices SOneTo(3)×SOneTo(3):
  0.368873   0.571601   0.732943
 -0.122903  -0.751632   0.648031
  0.921318  -0.329122  -0.207006

The direction obtained from the Minimum Variance Analysis (MVA) in our study closely aligns with the direction reported in the literature.

The (x, y, z) GSE components of the L, M, and N axes are L = (0.3665, –0.1201, 0.9226) GSE, M = (0.5694, –0.7553, –0.3245) GSE, and N = (0.7358, 0.6443, –0.2084) GSE

Plasma data

data_rate = "brst"
fpi_des_ds = FPIDataSet(; probe, data_rate, data_type="des")
fpi_dis_ds = FPIDataSet(; probe, data_rate, data_type="dis")

des_data = get_data(NamedTuple, fpi_des_ds, trange) |> dimarrayify
dis_data = get_data(NamedTuple, fpi_dis_ds, trange) |> dimarrayify

des_n = modify_meta(des_data.numberdensity[:, 1], label="Ne")
dis_bulkv_lmn = rotate(dis_data.bulkv_gse, rotMat) |> set_coord("LMN")
des_bulkv_lmn = rotate(des_data.bulkv_gse, rotMat) |> set_coord("LMN")
dis_n = modify_meta(dis_data.numberdensity[:, 1], label="Ni")
733-element DimArray{Unitful.Quantity{Float32, 𝐋^-3, Unitful.FreeUnits{(cm^-3,), 𝐋^-3, nothing}}, 1} mms2_dis_numberdensity_brst
├──────────────────────────────────────────────────────────────────────── dims ┤
  Ti Sampled{Dates.DateTime} [Dates.DateTime("2015-10-16T13:05:35.144"), …, Dates.DateTime("2015-10-16T13:07:24.946")] ForwardOrdered Irregular Points
├──────────────────────────────────────────────────────────────────── metadata ┤
  Dict{Any, Any} with 16 entries:
  "SCALETYP"        => "linear"
  "FILLVAL"         => Any[-1.0e31]
  "DEPEND_0"        => "Epoch"
  "FIELDNAM"        => "MMS2 FPI/DIS number density"
  "SI_CONVERSION"   => "1e6>m^-3"
  "VALIDMAX"        => Any[100000.0]
  :label            => "Ni"
  "DELTA_MINUS_VAR" => "mms2_dis_numberdensity_err_brst"
  "FORMAT"          => "E12.2"
  "VAR_TYPE"        => "data"
  "CATDESC"         => "MMS2 FPI/DIS ion number density during this burst"
  "LABLAXIS"        => "N"
  "DELTA_PLUS_VAR"  => "mms2_dis_numberdensity_err_brst"
  "DISPLAY_TYPE"    => "time_series"
  "VALIDMIN"        => Any[0.001]
  "UNITS"           => "cm^-3"
└──────────────────────────────────────────────────────────────────────────────┘
 ⋮      ⋱  

Magnitudes of electron and ion convection velocities

Vi_perp_mag, Ve_perp_mag = let Vi = dis_data.bulkv_gse, Ve = des_data.bulkv_gse, B = B_gse
    tr = timerange(B)
    Vi_clip = tview(Vi, tr)
    Ve_clip = tview(Ve, tr)

    B_int_Vi = tinterp(B, Vi_clip)[:, 1:3]
    B_int_Ve = tinterp(B, Ve_clip)[:, 1:3]
    Vi_perp_mag = toproj(Vi_clip, B_int_Vi) |> tnorm
    Ve_perp_mag = toproj(Ve_clip, B_int_Ve) |> tnorm

    modify_meta(Vi_perp_mag, label="Viper_t"),
    modify_meta(Ve_perp_mag, label="Veper_t")
end
(Quantity{Float32, 𝐋 𝐓^-1, Unitful.FreeUnits{(km, s^-1), 𝐋 𝐓^-1, nothing}}[16.370153f0 km s^-1, 39.0851f0 km s^-1, 18.778028f0 km s^-1, 36.386776f0 km s^-1, 28.834644f0 km s^-1, 31.232695f0 km s^-1, 21.860323f0 km s^-1, 25.041765f0 km s^-1, 19.903934f0 km s^-1, 3.4643996f0 km s^-1  …  54.849365f0 km s^-1, 52.696617f0 km s^-1, 43.563725f0 km s^-1, 37.47798f0 km s^-1, 32.163044f0 km s^-1, 29.341612f0 km s^-1, 22.57187f0 km s^-1, 23.173405f0 km s^-1, 26.802206f0 km s^-1, 33.881535f0 km s^-1], Quantity{Float32, 𝐋 𝐓^-1, Unitful.FreeUnits{(km, s^-1), 𝐋 𝐓^-1, nothing}}[179.42719f0 km s^-1, 65.01958f0 km s^-1, 95.88342f0 km s^-1, 80.233765f0 km s^-1, 31.386799f0 km s^-1, 176.10384f0 km s^-1, 21.706089f0 km s^-1, 53.19559f0 km s^-1, 126.49095f0 km s^-1, 102.78032f0 km s^-1  …  42.938835f0 km s^-1, 72.569305f0 km s^-1, 39.36957f0 km s^-1, 40.91847f0 km s^-1, 34.356808f0 km s^-1, 35.267933f0 km s^-1, 31.92348f0 km s^-1, 45.79629f0 km s^-1, 69.9509f0 km s^-1, 93.33416f0 km s^-1])

Current density

J = begin
    Vi = dis_bulkv_lmn
    Ve = des_bulkv_lmn
    n = dis_n
    Ve_clip = tclip(Ve, timerange(Vi))
    Vi_interp = tinterp(Vi, Ve_clip)
    n_interp = tinterp(n, Ve_clip)
    J = mapslices(Vi_interp - Ve_clip; dims=Ti) do V_diff_i
        @. Unitful.q * V_diff_i * n_interp |> u"μA*m^-2"
    end
    modify_meta(J, labels=["j_L", "j_M", "j_N"])
end
3661×3 DimArray{Unitful.Quantity{Float64, 𝐈 𝐋^-2, Unitful.FreeUnits{(μA, m^-2), 𝐈 𝐋^-2, nothing}}, 2}
├──────────────────────────────────────────────────────────────────────── dims ┤
  Ti Sampled{Dates.DateTime} [Dates.DateTime("2015-10-16T13:05:35.144"), …, Dates.DateTime("2015-10-16T13:07:24.946")] ForwardOrdered Irregular Points,
  AnonDim Sampled{Int64} 1:3 ForwardOrdered Regular Points
├──────────────────────────────────────────────────────────────────── metadata ┤
  Dict{Any, Any} with 19 entries:
  "UNITS"             => "km/s"
  "SCALETYP"          => "linear"
  "FILLVAL"           => Any[-1.0e31]
  "DEPEND_0"          => "Epoch"
  "FIELDNAM"          => "MMS2 FPI/DIS LMN bulk v"
  "LABL_PTR_1"        => ["Vx_LMN", "Vy_LMN", "Vz_LMN"]
  "SI_CONVERSION"     => "1.0e3>m s^-1"
  "VALIDMAX"          => Any[110000.0]
  "TENSOR_ORDER"      => Any[1]
  "COORDINATE_SYSTEM" => "LMN"
  "DELTA_MINUS_VAR"   => "mms2_dis_bulkv_err_brst"
  "FORMAT"            => "E12.2"
  "VAR_TYPE"          => "data"
  "CATDESC"           => "MMS2 FPI/DIS ion bulk-velocity LMN vector during this…
  :labels             => ["j_L", "j_M", "j_N"]
  "DELTA_PLUS_VAR"    => "mms2_dis_bulkv_err_brst"
  "DISPLAY_TYPE"      => "time_series"
  "VALIDMIN"          => Any[-110000.0]
  "REPRESENTATION_1"  => "mms2_dis_cartrep_brst"
└──────────────────────────────────────────────────────────────────────────────┘
 ⋮      ⋱  

Para, perp and anti-parallel spectra of electrons

des_energyspectr_tvars = [
    "cda/MMS2_FPI_BRST_L2_DES-MOMS/mms$(probe)_des_energyspectr_par_$(data_rate)",
    "cda/MMS2_FPI_BRST_L2_DES-MOMS/mms$(probe)_des_energyspectr_perp_$(data_rate)",
    "cda/MMS2_FPI_BRST_L2_DES-MOMS/mms$(probe)_des_energyspectr_anti_$(data_rate)"
]

des_energyspectr = get_data(des_energyspectr_tvars, trange) .|> DimArray
replace!.(des_energyspectr, 0 => NaN)
des_energyspectr_ratio = des_energyspectr[1] ./ des_energyspectr[2]
des_energyspectr = modify_meta.(des_energyspectr, colorrange=(1e5, 1e9))
des_energyspectr_ratio = modify_meta(des_energyspectr_ratio,
    colorrange=(1e-1, 1e1), title="Ratio (Para/Perp)"
)
[ Info: Cannot parse mms2_des_energyspectr_par_brst unit 
[ Info: Cannot parse mms2_des_energyspectr_perp_brst unit 
[ Info: Cannot parse mms2_des_energyspectr_anti_brst unit 
3666×32 DimArray{Float32, 2}
├──────────────────────────────┴───────────────────────────────────────── dims ┐
  Ti Sampled{Dates.DateTime} [Dates.DateTime("2015-10-16T13:05:35.024"), …, Dates.DateTime("2015-10-16T13:07:24.976")] ForwardOrdered Irregular Points,
  Y Sampled{Int64} 1:32 ForwardOrdered Regular Points
├──────────────────────────────────────────────────────────────────── metadata ┤
  Dict{Any, Any} with 20 entries:
  "VAR_NOTES"    => "Counts, summed within 30 degrees parallel bentPipe magneti…
  "SCALETYP"     => "log"
  "FILLVAL"      => Any[-1.0e31]
  "DEPEND_0"     => "Epoch"
  :ymeta         => PyDict{Any, Any}("CATDESC"=>"MMS FPI/DES burst sky-map pari…
  "FIELDNAM"     => "MMS2 FPI/DES energySpectr_par"
  "VALIDMAX"     => Any[1.0e30]
  :colorrange    => (0.1, 10.0)
  "DEPEND_1"     => "mms2_des_energy_brst"
  :yscale        => "log"
  "FORMAT"       => "E12.2"
  :ylabel        => "energy"
  "VAR_TYPE"     => "data"
  "CATDESC"      => "MMS2 FPI/DES electron energy parallel spectrum 30 degrees …
  "LABLAXIS"     => "FPI1/DES EnSpectr, Parallel"
  :title         => "Ratio (Para/Perp)"
  "DISPLAY_TYPE" => "spectrogram"
  "VALIDMIN"     => Any[-1.0e30]
  :yunit         => "eV"
  "y"            => mms2_des_energy_brst [ Units: eV, Size: (3666, 32)]
└──────────────────────────────────────────────────────────────────────────────┘
 ⋮      ⋱  

Summary data plot

In 1 paragraph (10 lines) explain what each panel represents for this reconnection interval.

From top to bottom, the panels represent the following: 1. the magnetic field vectors in LMN coordinate system; 2. energy times spectrogram of ion energy flux; 3. energy times spectrogram of electron energy flux; 4. plasma density (ion and electron); 5. ion flow velocity vectors in LMN coordinate system; 6. magnitudes of electron and ion convection velocities; 7. current density; 8. electron parallel and perpendicular temperatures; 9. electric field vectors in LMN coordinate system; 10-12. electron spectrograms (parallel, perpendicular, anti-parallel); 13. electron spectrogram ratio (para/perp).

tvars2plot = (
    B_LMN,
    [dis_n, des_n],
    modify_meta(dis_data.energyspectr_omni, colorrange=(1e3, 1e9)),
    modify_meta(des_data.energyspectr_omni, colorrange=(1e3, 1e9)),
    dis_bulkv_lmn,
    [Vi_perp_mag, Ve_perp_mag],
    E_LMN,
    des_energyspectr...,
    des_energyspectr_ratio
)
faxes = tplot(tvars2plot)
tlines!(faxes, "2015-10-16T13:05:52")
tlines!(faxes, ["2015-10-16T13:05:44", "2015-10-16T13:06:00"])
tlines!(faxes, ["2015-10-16T13:06:46", "2015-10-16T13:07:12"])
ylims!.(faxes.axes[end-3:end], 0.9e1, 1.2e3)
faxes
Warning: y values are not constant along time
@ SPEDAS ~/.julia/packages/SPEDAS/TYwzF/src/meta.jl:46
Warning: y values are not constant along time
@ SPEDAS ~/.julia/packages/SPEDAS/TYwzF/src/meta.jl:46
Warning: y values are not constant along time
@ SPEDAS ~/.julia/packages/SPEDAS/TYwzF/src/meta.jl:46
Warning: y values are not constant along time
@ SPEDAS ~/.julia/packages/SPEDAS/TYwzF/src/meta.jl:46
Warning: y values are not constant along time
@ SPEDAS ~/.julia/packages/SPEDAS/TYwzF/src/meta.jl:46
Warning: y values are not constant along time
@ SPEDAS ~/.julia/packages/SPEDAS/TYwzF/src/meta.jl:46
Warning: y values are not constant along time
@ SPEDAS ~/.julia/packages/SPEDAS/TYwzF/src/meta.jl:46
Warning: y values are not constant along time
@ SPEDAS ~/.julia/packages/SPEDAS/TYwzF/src/meta.jl:46
Warning: y values are not constant along time
@ SPEDAS ~/.julia/packages/SPEDAS/TYwzF/src/meta.jl:46
Warning: y values are not constant along time
@ SPEDAS ~/.julia/packages/SPEDAS/TYwzF/src/meta.jl:46
Warning: y values are not constant along time
@ SPEDAS ~/.julia/packages/SPEDAS/TYwzF/src/meta.jl:46
Warning: y values are not constant along time
@ SPEDAS ~/.julia/packages/SPEDAS/TYwzF/src/meta.jl:46
Back to top