pyrfu.pyrf.ebsp module#

pyrfu.pyrf.ebsp.ebsp(e_xyz, db_xyz, b_xyz, b_bgd, xyz, freq_int, **kwargs)[source]#

Calculates wavelet spectra of E&B and Poynting flux using wavelets (Morlet wavelet). Also computes polarization parameters of B using SVD [7]. SVD is performed on spectral matrices computed from the time series of B using wavelets and then averaged over a number of wave periods.

Parameters:
  • e_xyz (xarray.DataArray) – Time series of the wave electric field.

  • db_xyz (xarray.DataArray) – Time series of the wave magnetic field.

  • b_xyz (xarray.DataArray) – Time series of the high resolution background magnetic field used for E.B=0.

  • b_bgd (xarray.DataArray) – Time series of the background magnetic field used for field aligned coordinates.

  • xyz (xarray.DataArray) – Time series of the position time series of spacecraft used for field aligned coordinates.

  • freq_int (str or array_like) –

    Frequency interval :
    • ”pc12” : [0.1, 5.0]

    • ”pc35” : [2e-3, 0.1]

    • [fmin, fmax] : arbitrary interval [fmin,fmax]

  • polarization (bool) – Computes polarization parameters. Default False.

  • no_resample (bool) – No resampling, E and delta_b are given at the same time line. Default False.

  • fac (bool) – Uses FAC coordinate system (defined by b0 and optionally xyz), otherwise no coordinate system transformation is performed. Default True.

  • de_dot_b0 (bool) – Computes dEz from delta_b dot B = 0, uses full_b. Default False.

  • full_b_db (bool) – delta_b contains DC field. Default False.

  • nav (int) – Number of wave periods to average Default 8.

  • fac_matrix (numpy.ndarray) – Specify rotation matrix to FAC system Default None.

  • m_width_coeff (int or float) – Specify coefficient to multiple Morlet wavelet width by. Default 1.

Returns:

res

Dataset with :
  • txarray.DataArray

    Time.

  • fxarray.DataArray

    Frequencies.

  • bb_xxyyzzssxarray.DataArray
    delta_b power spectrum with :
    • […,0] : x

    • […,1] : y

    • […,2] : z

    • […,3] : sum

  • ee_xxyyzzssxarray.DataArray
    E power spectrum with :
    • […,0] : x

    • […,1] : y

    • […,2] : z

    • […,3] : sum

  • ee_ssxarray.DataArray

    E power spectrum (xx+yy spacecraft coordinates, e.g. ISR2).

  • pf_xyzxarray.DataArray

    Poynting flux (xyz).

  • pf_rtpxarray.DataArray

    Poynting flux (r, theta, phi) [angles in degrees].

  • dopxarray.DataArray

    3D degree of polarization.

  • dop2dxarray.DataArray

    2D degree of polarization in the polarization plane.

  • planarityxarray.DataArray

    Planarity of polarization.

  • ellipticityxarray.DataArray

    Ellipticity of polarization ellipse.

  • kxarray.DataArray

    k-vector (theta, phi FAC) [angles in degrees].

Return type:

xarray.Dataset

See also

pyrfu.plot.pl_ebsp

to fill.

pyrfu.pyrf.convert_fac

Transforms to a field-aligned coordinate.

Notes

This software was developed as part of the MAARBLE (Monitoring, Analyzing and Assessing Radiation Belt Energization and Loss) collaborative research project which has received funding from the European Community’s Seventh Framework Programme (FP7-SPACE-2011-1) under grant agreement n. 284520.

References

[7]

Santolík, O., Parrot. M., and Lefeuvre. F. (2003) Singular value decomposition methods for wave propagation analysis,Radio Sci., 38(1), 1010, doi : https://doi.org/10.1029/2000RS002523 .

Examples

>>> from pyrfu import mms, pyrf
>>> # Time interval
>>> tint_brst = ["2015-10-30T05:15:42.000", "2015-10-30T05:15:54.000"]
>>> # Spacecraft index
>>> mms_id = 3
>>> # Load spacecraft position
>>> tint_long = pyrf.extend_tint(tint_brst, [-100, 100])
>>> r_xyz = mms.get_data("R_gse", tint_long, mms_id)
>>> # Load background magnetic field, electric field and magnetic field
fluctuations
>>> b_xyz = mms.get_data("B_gse_fgm_brst_l2", tint_brst, mms_id)
>>> e_xyz = mms.get_data("E_gse_edp_brst_l2", tint_brst, mms_id)
>>> b_scm = mms.get_data("B_gse_scm_brst_l2", tint_brst, mms_id)
>>> # Polarization analysis
>>> options = dict(polarization=True, fac=True)
>>> polarization = pyrf.ebsp(e_xyz, b_scm, b_xyz, b_xyz, r_xyz,
>>>                          freq_int=[10, 4000], **options)