SwarmFACE package
Subpackages
- SwarmFACE.plot_save package
- Submodules
- SwarmFACE.plot_save.dual_sat_BI module
- SwarmFACE.plot_save.dual_sat_LS module
- SwarmFACE.plot_save.dual_sat_SVD module
- SwarmFACE.plot_save.saveplot_qi module
- SwarmFACE.plot_save.single_sat module
- SwarmFACE.plot_save.single_sat_MVA module
- SwarmFACE.plot_save.three_sat module
- SwarmFACE.plot_save.three_sat_conj module
- Module contents
Submodules
SwarmFACE.esaL2 module
- SwarmFACE.esaL2.dual(dtime_beg, dtime_end)[source]
Retrieve the Level-2 dual-satellite FAC density from the ESA database
- Parameters
dtime_beg (str) – start time in ISO format ‘YYYY-MM-DDThh:mm:ss’
dtime_end (str) – end time in ISO format
- Returns
FAC_L2 – the Level-2 dual-satellite FAC and IRC densities
- Return type
DataFrame
- SwarmFACE.esaL2.single(dtime_beg, dtime_end, sat)[source]
Retrieve the Level-2 single-satellite FAC density from the ESA database
- Parameters
dtime_beg (str) – start time in ISO format ‘YYYY-MM-DDThh:mm:ss’
dtime_end (str) – end time in ISO format
sat ([str]) – satellite, e.g. [‘A’]
- Returns
FAC_L2 – the Level-2 single-satellite FAC and IRC densities
- Return type
DataFrame
SwarmFACE.fac module
- SwarmFACE.fac.bi_dualJfac(dt, R, B, dB, dt_along=5, er_db=0.5, angTHR=30.0, errTHR=0.1, use_filter=True, saveconf=False)[source]
Estimate the FAC density with dual-satellite Boundary Integral method
- Parameters
dt (pandas.Index) – time stamps
R (numpy.array) – satellites vector position in GEO frame
B (numpy.array) – magnetic field vector in GEO frame
dB (numpy.array) – magnetic perturbation vector in GEO frame
dt_along (integer) – along track separation in sec.
er_db (float) – error in magnetic field measurements
angTHR (float) – minimum accepted angle between the magnetic field vector and the quad plane
errTHR (float) – accepted error for normal current density
use_filter (boolean) – ‘True’ for data filtering
saveconf (boolean) – ‘True’ to add the quad’s parameters in the results
- Returns
j_df – results, including FAC and normal current densities
- Return type
DataFrame
- SwarmFACE.fac.calcBI(R4, dB4)[source]
Compute curlB using the discrete BI integral along a four-point configuration.
- Parameters
R4 (numpy.array) – position vectors of the apexes in format [index, apex, comp]
dB4 (numpy.array) – magnetic perturbation vectors at the apexes in format [index, apex, comp]
- Returns
curl of magnetic field perturbation
- Return type
numpy.array
- SwarmFACE.fac.ffunct(tau3d, tauast=0.13, taunul=0.07, intpol='Linear')[source]
Compute coefficients for a smooth transition between 1D and 2D gradient estimators in the SVD method.
- Parameters
tau3d (numpy.array) – S eigenvalues ratio
tauast (float) – value between [0, 1]. tauast >= taunul. The interval [taunul, tauast] is used to match the 1D and 2D gradient estimators
taunul (float) – value between [0, 1]. Below this threshold value the quad is considered degenerated
intpol (string) – interpolation method used for matching the 1D and 2D gradient estimators. Should be ‘Linear’, ‘Cubic’, or None
- Returns
f – coefficient to smooth transition between 1D and 2D gradient estimators
- Return type
numpy.array
- SwarmFACE.fac.ls_dualJfac(dt, R, B, dB, dt_along=5, er_db=0.5, angTHR=30.0, errTHR=0.1, use_filter=True, saveconf=False)[source]
Estimate the FAC density with dual-satellite Least-Squares method
- Parameters
dt (pandas.Index) – time stamps
R (numpy.array) – satellites vector position in GEO frame
B (numpy.array) – magnetic field vector in GEO frame
dB (numpy.array) – magnetic perturbation vector in GEO frame
dt_along (integer) – along track separation in sec.
er_db (float) – error in magnetic field measurements
angTHR (float) – minimum accepted angle between the magnetic field vector and the quad plane
errTHR (float) – accepted error for normal current density
use_filter (boolean) – ‘True’ for data filtering
saveconf (boolean) – ‘True’ to add the quad’s parameters in the results
- Returns
j_df – results, including FAC and normal current densities
- Return type
DataFrame
- SwarmFACE.fac.qmatrix_intpol(R, tauast=0.13, taunul=0.07, intpol='Linear')[source]
Compute the SVD reciprocal matrix of a four-point planar configuration. A smooth transition between 1D and 2D gradient estimators is allowed.
- Parameters
R (numpy.array) – position vectors of the apexes in format [index, apex, comp]
tauast (float) – value between [0, 1]. tauast >= taunul. The interval [taunul, tauast] is used to match the 1D and 2D gradient estimators
taunul (float) – value between [0, 1]. Below this threshold value the quad is considered degenerated
intpol (string) – interpolation method used for matching the 1D and 2D gradient estimators. Should be ‘Linear’, ‘Cubic’, or None
- Returns
Q (numpy.array) – reciprocal matrix
AD (numpy_array) – quad degeneracy
Vt (numpy.array) – orthogonal matrix from SVD decomposition of R = USVt
S (numpy.array) – diagonal matrix from SVD decomposition of R = USVt
- SwarmFACE.fac.recivec2s(R4s)[source]
Compute the reciprocal vectors of a four-point planar configuration
- Parameters
R4s (numpy.array) – position vectors of the apexes in format [index, apex, comp]
- Returns
Q4s (numpy.array) – canonical base (reciprocal) vectors
Rc (numpy_array) – position of the mesocenter
nuvec (numpy.array) – normal to the spacecraft plane
condnum (numpy.array) – condition number
nonplanar (numpy.array) – nonplanarity
tracer (numpy.array) – trace of the position tensor
traceq (numpy.array) – trace of the reciprocal tensor
- SwarmFACE.fac.recivec3s(R)[source]
Compute the reciprocal vectors for a three-satellites configuration
- Parameters
R (numpy.array) – satellites vector position in format [index, sc, comp]
- Returns
Rc (numpy.array) – position of the mesocenter
Rmeso (numpy.array) – satellites position in the mesocentric frame
nuvec (numpy.array) – normal to the spacecraft plane
Q3S (numpy.array) – generalized reciprocal vectors
Qtens (numpy.array) – reciprocal tensor
Rtens (numpy.array) – position tensor
- SwarmFACE.fac.singleJfac(t, R, B, dB, alpha=None, N2d=None, N3d=None, tincl=None, res='LR', er_db=0.5, angTHR=30.0, use_filter=True)[source]
Estimate the FAC density with single-satellite method
- Parameters
t (pandas.Index) – time stamps
R (numpy.array) – satellite vector position in GEO frame
B (numpy.array) – magnetic field vector in GEO frame
dB (numpy.array) – magnetic perturbation vector in GEO frame
alpha (float) – angle in the tangential plane between the (projection of) current sheet normal and the satellite velocity
N3d ([float, float, float]) – current sheet normal vector in GEO frame
N2d ([float, float, float]) – projection of the current sheet normal on the tangential plane
tincl ([datetime, datetime]) – time interval when the information on current sheet inclination is valid
res (str) – data resolution, ‘LR’ or ‘HR’
er_db (float) – error in magnetic field measurements
angTHR (float) – minimum accepted angle between the magnetic field vector and the tangential plane
use_filter (boolean) – ‘True’ for data filtering
- Returns
j_df – results, including FAC and IRC densities
- Return type
DataFrame
- SwarmFACE.fac.svd_dualJfac(dt, R, B, dB, dt_along=5, er_db=0.5, tauast=0.13, taunul=0.07, intpol='Linear', angTHR=30.0, use_filter=True, saveconf=False)[source]
Estimate the FAC density with dual-satellite Singular Values Decomposition method
- Parameters
dt (pandas.Index) – time stamps
R (numpy.array) – satellites vector position in GEO frame
B (numpy.array) – magnetic field vector in GEO frame
dB (numpy.array) – magnetic perturbation vector in GEO frame
dt_along (integer) – along track separation in sec.
er_db (float) – error in magnetic field measurements
tauast (float) – value between [0, 1]. tauast >= taunul. The interval [taunul, tauast] is used to match the 1D and 2D gradient estimators
taunul (float) – value between [0, 1]. Below this threshold value the quad is considered degenerated
intpol (string) – interpolation method used for matching the 1D and 2D gradient estimators. Should be ‘Linear’, ‘Cubic’, or None
angTHR (float) – minimum accepted angle between the magnetic field vector and the quad plane
use_filter (boolean) – ‘True’ for data filtering
saveconf (boolean) – ‘True’ to add the quad’s parameters in the results
- Returns
j_df – results, including FAC and normal current densities
- Return type
DataFrame
- SwarmFACE.fac.threeJfac(dt, R, B, dB, er_db=0.5, angTHR=30.0, use_filter=True)[source]
Estimate the FAC density with three-satellite method
- Parameters
dt (pandas.Index) – time stamps
R (numpy.array) – satellites vector position in GEO frame
B (numpy.array) – magnetic field vector in GEO frame
dB (numpy.array) – magnetic perturbation vector in GEO frame
er_db (float) – error in magnetic field measurements
angTHR (float) – minimum accepted angle between the magnetic field vector and the spacecraft plane
use_filter (boolean) – ‘True’ for data filtering
- Returns
j_df – results, including FAC and normal current densities
- Return type
DataFrame
SwarmFACE.fac_qi module
- SwarmFACE.fac_qi.fac_qi(dtime_beg, dtime_end, swB=False, jTHR=0.05, saveplot=False)[source]
High-level routine to automatically estimate the quality indices of FAC structures
- Parameters
dtime_beg (str) – start time in ISO format ‘YYYY-MM-DDThh:mm:ss’
dtime_end (str) – end time in ISO format
swB (boolean) – ‘True’ to the MVA based quality indices for Swarm B
jTHR (float) – threshold to neglect small FAC densities
saveplot (boolean) – ‘True’ for plotting the results
- Returns
input_df (list) – list of DataFrames (one per satellite per quarter-orbit section) with input data
RBdBAng_df (list) – list of DataFrrames (one per satellite per quarter-orbit section) with intermediate variables
fac_df (list) – list of DataFrrames (one per satellite per quarter-orbit section) with estimated single=satellite FAC data
qimva_df (DataFrame) – MVA results, FAC planarity, and FAC inclination
qicc_df (DataFrame) – results from the correlation analysis between the magnetic field perturbations at the lower Swarm satellites
param (dict) – parameters used in the analysis
- SwarmFACE.fac_qi.qorbsCC(sats, qimva_df, qorbs_dBmva)[source]
Perform correlation analysis between the magnetic field perturbations at the lower Swarm satellite
- Parameters
sats (list) – satellites entering the analysis, e.g. [‘A’, ‘C’] or [‘A’, ‘B’, ‘C’]
qimva_df (DataFrame) – MVA results, FAC planarity, and FAC inclination
qorbs_dBmva (list) – list of DataFrrames (one per satellite per quarter-orbit section) with magnetic field perturbation in MVA frame
- Returns
qicc_df (DataFrame) – results from the correlation analysis between the magnetic field perturbations at the lower Swarm satellites
Bcc_df (list) – list of DataFrrames (one per quarter-orbit section) with magnetic field perturbation of the reference satellite along the maximum magnetic variance direction
- SwarmFACE.fac_qi.qorbsMVA(sats, orbs, qorbs_dB, qorbs_fac)[source]
Automatically estimate the MVA related quality indices of FAC structures
- Parameters
sats (list) – satellites entering the analysis, e.g. [‘A’, ‘C’] or [‘A’, ‘B’, ‘C’]
orbs (array) – array of integers of size [nr. sat, 2], indicate for each satellite the start and end orbit numbers
qorbs_dB (list) – list of DataFrrames (one per satellite per quarter-orbit section) with magnetic field perturbation in GEO frame
fac_df (list) – list of DataFrrames (one per satellite per quarter-orbit section) with estimated single=satellite FAC data
- Returns
qimva_df (DataFrame) – MVA results, FAC planarity, and FAC inclination
qorbs_dBmva (list) – list of DataFrrames (one per satellite per quarter-orbit section) with magnetic field perturbation in MVA frame
SwarmFACE.find_3sat_conj module
- SwarmFACE.find_3sat_conj.find_3sat_conj(dtime_beg, dtime_end, delT=120, delN=300, delE=1200, jTHR=0.05, saveplot=True)[source]
High-level routine to find Swarm conjunctions above the auroral oval
- Parameters
dtime_beg (str) – start time in ISO format ‘YYYY-MM-DDThh:mm:ss’
dtime_end (str) – end time in ISO format
delT (float) – time window (in sec.) to find conjunctions
delN (float) – spatial range along North-South (in km) to find conjunctions
delE (float) – spatial range along East-West (in km) to find conjunctions
jTHR (float) – threshold to neglect small FAC densities
saveplot (boolean) – ‘True’ for plotting the results
- Returns
conj_df (DataFrame) – details on the identified conjunctions
param (dict) – parameters used in the analysis
SwarmFACE.j1sat module
- SwarmFACE.j1sat.j1sat(dtime_beg, dtime_end, sat, res='LR', use_filter=True, alpha=None, N3d=None, N2d=None, tincl=None, er_db=0.5, angTHR=30.0, savedata=True, saveplot=True)[source]
High-level routine to estimate the FAC density with single-satellite method
- Parameters
dtime_beg (str) – start time in ISO format ‘YYYY-MM-DDThh:mm:ss’
dtime_end (str) – end time in ISO format
sat ([str]) – satellite, e.g. [‘A’]
res (str) – data resolution, ‘LR’ or ‘HR’
use_filter (boolean) – ‘True’ for data filtering
alpha (float) – angle in the tangential plane between the (projection of) current sheet normal and the satellite velocity
N3d ([float, float, float]) – current sheet normal vector in GEO frame
N2d ([float, float, float]) – projection of the current sheet normal on the tangential plane
tincl ([ISO time, ISO time]) – time interval when the information on current sheet inclination is valid
er_db (float) – error in magnetic field measurements
angTHR (float) – minimum accepted angle between the magnetic field vector and the tangential plane
savedata (boolean) – ‘True’ for saving the results in an ASCII file
saveplot (boolean) – ‘True’ for plotting the results
- Returns
j_df (DataFrame) – results
input_df (DataFrame) – input data
param (dict) – parameters used in the analysis
SwarmFACE.j2satBI module
- SwarmFACE.j2satBI.j2satBI(dtime_beg, dtime_end, sats, tshift=None, dt_along=5, use_filter=True, er_db=0.5, angTHR=30.0, errTHR=0.1, saveconf=False, savedata=True, saveplot=True)[source]
High-level routine to estimate the FAC density with dual-satellite Boundary Integral method
- Parameters
dtime_beg (str) – start time in ISO format ‘YYYY-MM-DDThh:mm:ss’
dtime_end (str) – end time in ISO format
sats ([str, str]) – satellite pair, e.g. [‘A’, ‘C’]
tshift ([float, float]) – array of time shifts (in seconds) to be introduced in satellite data in order to achieve the desired quad configuration
dt_along (int) – quad’s length in the along track direction (in seconds of satellite travel distance)
use_filter (boolean) – ‘True’ for data filtering
er_db (float) – error in magnetic field measurements
angTHR (float) – minimum accepted angle between the magnetic field vector and the quad’s plane
errTHR (float) – accepted error for the current density along the normal direction
saveconf (boolean) – ‘True’ for adding the quad’s parameters in the results
savedata (boolean) – ‘True’ for saving the results in an ASCII file
saveplot (boolean) – ‘True’ for plotting the results
- Returns
j_df (DataFrame) – results
input_df (DataFrame) – input data
param (dict) – parameters used in the analysis
SwarmFACE.j2satLS module
- SwarmFACE.j2satLS.j2satLS(dtime_beg, dtime_end, sats, tshift=None, dt_along=5, use_filter=True, er_db=0.5, angTHR=30.0, errTHR=0.1, saveconf=False, savedata=True, saveplot=True)[source]
High-level routine to estimate the FAC density with dual-satellite Least Squares method
- Parameters
dtime_beg (str) – start time in ISO format ‘YYYY-MM-DDThh:mm:ss’
dtime_end (str) – end time in ISO format
sats ([str, str]) – satellite pair, e.g. [‘A’, ‘C’]
tshift ([float, float]) – array of time shifts (in seconds) to be introduced in satellite data in order to achieve the desired quad configuration
dt_along (int) – quad’s length in the along track direction (in seconds of satellite travel distance)
use_filter (boolean) – ‘True’ for data filtering
er_db (float) – error in magnetic field measurements
angTHR (float) – minimum accepted angle between the magnetic field vector and the quad’s plane
errTHR (float) – accepted error for the current density along the normal direction
saveconf (boolean) – ‘True’ for adding the quad’s parameters in the results
savedata (boolean) – ‘True’ for saving the results in an ASCII file
saveplot (boolean) – ‘True’ for plotting the results
- Returns
j_df (DataFrame) – results
input_df (DataFrame) – input data
param (dict) – parameters used in the analysis
SwarmFACE.j2satSVD module
- SwarmFACE.j2satSVD.j2satSVD(dtime_beg, dtime_end, sats, tshift=None, dt_along=5, use_filter=True, er_db=0.5, tauast=0.13, taunul=0.07, intpol='Linear', angTHR=30.0, saveconf=False, savedata=True, saveplot=True)[source]
High-level routine to estimate the FAC density with dual-satellite Singular Value Decomposition method
- Parameters
dtime_beg (str) – start time in ISO format ‘YYYY-MM-DDThh:mm:ss’
dtime_end (str) – end time in ISO format
sats ([str, str]) – satellite pair, e.g. [‘A’, ‘C’]
tshift ([float, float]) – array of time shifts (in seconds) to be introduced in satellite data in order to achieve the desired quad configuration
dt_along (int) – quad’s length in the along track direction (in seconds of satellite travel distance)
use_filter (boolean) – ‘True’ for data filtering
er_db (float) – error in magnetic field measurements
angTHR (float) – minimum accepted angle between the magnetic field vector and the quad’s plane
tauast (float) – value between [0, 1]. tauast >= taunul. The interval [taunul, tauast] is used to match the 1D and 2D gradient estimators
taunul (float) – value between [0, 1]. Below this threshold value the quad is considered degenerated
intpol (str) – interpolation method used for matching the 1D and 2D gradient estimators. Should be ‘Linear’, ‘Cubic’, or None
saveconf (boolean) – ‘True’ for adding the quad’s parameters in the results
savedata (boolean) – ‘True’ for saving the results in an ASCII file
saveplot (boolean) – ‘True’ for plotting the results
- Returns
j_df (DataFrame) – results
input_df (DataFrame) – input data
param (dict) – parameters used in the analysis
SwarmFACE.j3sat module
- SwarmFACE.j3sat.j3sat(dtime_beg, dtime_end, tshift=[0, 0, 0], use_filter=True, er_db=0.5, angTHR=30.0, savedata=True, saveplot=True)[source]
High-level routine to estimate the FAC density with three-satellite method
- Parameters
dtime_beg (str) – start time in ISO format ‘YYYY-MM-DDThh:mm:ss’
dtime_end (str) – end time in ISO format
tshift ([float, float, float]) – array of time shifts (in seconds) to be introduced in satellite data in order to achieve a more favorable configuration
use_filter (boolean) – ‘True’ for data filtering
er_db (float) – error in magnetic field measurements
angTHR (float) – minimum accepted angle between the magnetic field vector and the spacecraft plane
savedata (boolean) – ‘True’ for saving the results in an ASCII file
saveplot (boolean) – ‘True’ for plotting the results
- Returns
j_df (DataFrame) – results
input_df (DataFrame) – input data
param (dict) – parameters used in the analysis
SwarmFACE.mva1sat module
- SwarmFACE.mva1sat.get_data_mva1sat(dtime_beg, dtime_end, sat, use_filter=True)[source]
Prepare and plot data on the screen to interactively select the MVA interval
- Parameters
dtime_beg (str) – start time in ISO format ‘YYYY-MM-DDThh:mm:ss’
dtime_end (str) – end time in ISO format
sat ([str]) – satellite, e.g. [‘A’]
use_filter (boolean) – ‘True’ for data filtering
- Returns
j_df (DataFrame) – FAC density from single-satellite method
input_df (DataFrame) – input data (includes magnetic field perturbation)
param (dict) – parameters used in the analysis
span_sel (dict) – new MVA start and stop time
span (reference) – reference to SpanSelector to prevent it from being garbage collected
- SwarmFACE.mva1sat.perform_mva1sat(j_df, input_df, param, span_sel, savedata=True, saveplot=True)[source]
Perform the MVA analysis on the data provided by get_data_mva1sat
- Parameters
j_df (DataFrame) – FAC density from single-satellite method
input_df (DataFrame) – input data (includes magnetic field perturbation)
param (dict) – parameters used in the analysis
span_sel (dict) – new MVA start and stop time
savedata (boolean) – ‘True’ for saving the results in an ASCII file
saveplot (boolean) – ‘True’ for plotting the results
- Returns
jcorr_df (DataFrame) – corrected (for current sheet inclination) FAC density
dBmva_df (DataFrame) – magnetic field perturbation in MVA frame
mva_df (DataFrame) – MVA results, FAC planarity and inclination
param (dict) – parameters used in the analysis
SwarmFACE.utils module
- SwarmFACE.utils.GapsAsNaN(df_ini, ind_gaps)[source]
Given a DataFrame object, replace with NaN the magnetic field data and model for certain location
- Parameters
df_ini (DataFrame) – containes columns with magetic field data and model
ing_gaps (numpy.array) – integer-location based indexing
- Returns
df_out (DataFrame) – new object with NaN
index – index values at replaced location
- SwarmFACE.utils.R_in_GEOC(Rsph)[source]
Returns Swarm position in cartesian GEO frame and the rotation matrix from NEC to GEO
- Parameters
Rsph (numpy.array) – satellite position in spherical GEO frame
- Returns
R (numpy.array) – satellite position in cartesian GEO frame
numpy.array – array of rotation matrices from NEC to GEO
- SwarmFACE.utils.SortVertices(R4s, dB4s)[source]
Sort the quad’s vertices in correct order; compute quad’s parameters and the trace of the BI and LS reciprocal tensors.
- Parameters
R4s (numpy.array) – position vectors of the apexes in format [index, apex, comp]
dB4s (numpy.array) – magnetic perturbation vectors at the apexes in format [index, apex, comp]
- Returns
R4sa, dB4sa (numpy.arrays) – position and magnetic perturbation vectors sorted in correct order
trBI (numpy.array) – trace of the reciprocal tensor in LS method
Rmeso (numpy.array) – position of the mesocenter
nuvec (numpy.array) – normal to the quad
satsort (numpy.array) – satellite indeces in the correct order
trLS (numpy.array) – trace of the reciprocal tensor in LS method
EL, EM, el, em (numpy.arrays) – parameters of the four-point configuration
- SwarmFACE.utils.find_ao_margins(df, fac_qnt='FAC_flt_sup', rez_qd=100)[source]
Estimate the times of auroral oval (AO) encounter (central point and edges) for a quarte-orbit section. Works with quasi-dipole latitude QDLat (not time) dependence.
- Parameters
df (DataFrame) – Swarm data for a quarter-orbit section, including FAC data
fac_qnt (string) – name of the column with FAC data to work with
rez_qd – resolution to interpolate the QDLat
- Returns
t64[0], t64[1], t64[2] (datetime-like) – AO start, stop, and central time
float, float – QDLat and QDLon of AO center
qd_trend, qd_sign (float) – trend and sign of QDLat in the analysed quarte-orbit section
ti_arr, qd_arr (numpy.array) – interpolated times and corresponding QDLat values
- SwarmFACE.utils.find_tshift2sat(dti, Rsph, sats)[source]
Compute orbital phase lag [in sec.] between two satellites
- Parameters
dti (Datetimeindex) – One second separated timestamps
Rshp (numpy.array) – satellites position in spherical GEO frame
sats ([str, str]) – satellite pair, e.g. [‘A’, ‘C’]
- Returns
tshift – optimal time shift in sec. to aligned the two sensors side by side
- Return type
[float, float]
- SwarmFACE.utils.get_db_data(sats, tbeg, tend, Bmodel, frame='GEOC')[source]
Return Swarm magnetic perturbation in a DataFrame structure
- Parameters
sats (str or list of str) – satellites
tbeg (datetime-like) – start and stop times
tend (datetime-like) – start and stop times
Bmodel (string) – magnetic field model used to compute the magnetic field perturbations (currently CHAOS)
frame (string) – reference frame for the output, i.e. ‘GEOC’ for GEO frame or None for NEC frame
- Returns
dB – Swarm magnetic field perturbation data
- Return type
DataFrame
- SwarmFACE.utils.mva(v, cdir=None)[source]
Return results of Minimum Variance Analysis on an array of 3D vector v. If cdir (3D vector) is provided, the analysis is performed in the plane perpendiculat to it
- Parameters
v (numpy.array) – 3D vectors
cdir (numpy.array) – direction to constrain MVA
- Return type
array with eigenvalues and eigenvectors from MVA
- SwarmFACE.utils.normvec(v)[source]
Return array of unit vectors
- Parameters
v (numpy.array) – vectors
- Return type
numpy.array of unit vectors
- SwarmFACE.utils.res_param(res)[source]
Return sampling step and sampling frequency acording to Swarm magnetic data resolution
- Parameters
res (string) – ‘HR’ or ‘LR’
- Return type
sampling step (string) and sampling frequency (int)
- SwarmFACE.utils.rotvecax(v, ax, ang)[source]
Rotates vector v by angle ang around a normal vector ax. Uses Rodrigues’ formula when v is normal to ax
- Parameters
v (numpy.array) – vector(s)
ax (numpy.array) – normal vector(s)
ang (numpy.array) – angle(s)
- Return type
vector or array of vectors
- SwarmFACE.utils.sign_ang(V, N, R)[source]
Return the signed angle between vectors V and N, perpendicular to R; positive sign corresponds to right hand rule along R
- Parameters
V (numpy.array) – vectors
N (numpy.array) – vectors
R (numpy.array) – vectors
- Return type
angle(s)
- SwarmFACE.utils.sign_dang(ang2, ang1)[source]
Return the difference between two angles in [-180, 180] deg. interval
- Parameters
ang1 (numpy.array) – angles in degrees
ang2 (numpy.array) – angles in degrees
- Return type
numpy.array of angles
- SwarmFACE.utils.split_into_sections(df, begend_arr)[source]
Splits the original DataFrame object df in more DataFrame objects according to time moments from begend_arr array.
- Parameters
df (DataFrame) – time indexed DataFrame object
begend_arr (list) – time moments
Returns – list of DataFrame objects