swarmpal.toolboxes.dsecs#
swarmpal.toolboxes.dsecs.processes#
- class swarmpal.toolboxes.dsecs.processes.Analysis(config: dict | None = None, active_tree: str = '/', inplace: bool = True)[source]#
Bases:
PalProcess
Run the DSECS analysis
- property active_tree: str#
Defines which branch of the datatree will be used
- property config: dict#
Dictionary that configures the process behaviour
- property process_name#
- class swarmpal.toolboxes.dsecs.processes.Preprocess(config: dict | None = None, active_tree: str = '/', inplace: bool = True)[source]#
Bases:
PalProcess
Prepare data for input to DSECS analysis
- property active_tree: str#
Defines which branch of the datatree will be used
- property config: dict#
Dictionary that configures the process behaviour
- property process_name#
swarmpal.toolboxes.dsecs.dsecs_algorithms#
Algorithms for low latitude spherical elementary current system analysis.
Adapted from MatLab code by Heikki Vanhamäki.
- swarmpal.toolboxes.dsecs.dsecs_algorithms.SECS_1D_CurlFree_magnetic(latB, latSECS, rb, rsecs, geometry)[source]#
Calculates the array that maps the CF SECS amplitudes to magnetic field.
- Parameters:
latB (ndarray) – Latitudes where the magnetic field is calculated.
latSECS (ndarray) – Latitudes of the SECS poles
rb (ndarray) – Array of radial distances of the points where the magnetic field is calculated, scalar or vector, [km]
rsecs (float) – Radius of the sphere where the calculation takes place, [km], scalar
geometry (string) – If ‘radial’, assume radial field lines. Else, assume dipolar field.
- Returns:
2D array that maps the SECS amplitudes to magnetic field measurements.
- Return type:
ndarray
- swarmpal.toolboxes.dsecs.dsecs_algorithms.SECS_1D_CurlFree_vector(latV, latSECS, ri)[source]#
Calculates the array relating SECS amplitudes to curl free current.
- Parameters:
latV (ndarray) – Array of latitudes where the current is calculated.
latSECS (ndarray) – Array of latitudes of the SECS poles.
ri (float) – Assumed radius of the spherical shell where the currents flow (km).
- Returns:
2D array to map the elementary current amplitudes to div free current.
- Return type:
ndarray
- swarmpal.toolboxes.dsecs.dsecs_algorithms.SECS_1D_DivFree_magnetic(latB, latSECS, rb, rsecs, M)[source]#
Calculates the matrix for magnetic field from 1D divergence free current.
- Parameters:
latB (ndarray) – Array of latitudes where the magnetic field is calculated [degree]
latSECS (ndarray) – Array of latitudes where 1D CF SECS are located [degree]
rb (ndarray) – Array of radial distances of the points where the magnetic field is calculated, scalar or vector, [km]
rsecs (float) – Radius of the sphere where the 1D CF SECS are located, scalar, [km]
M (int) – argest order of the Legendre polynomials
- Returns:
ndarray – 2D array to map the elementary current system amplitudes to radial magnetic field
ndarray – 2D array to map the elementary current system amplitudes to meridional magnetic field
- swarmpal.toolboxes.dsecs.dsecs_algorithms.SECS_1D_DivFree_vector(latV, latSECS, ri)[source]#
Calculates the matrix for 1D divergence free current density.
- Parameters:
latV (ndarray) – Latitudes of locations where vectors is calculated, in degrees
latSECS (ndarray) – Latitudes of the 1D SECS poles, in degrees
ri (int) – Assumed radius of the spherical shell where the currents flow (km).
- Returns:
2D array that maps the 1D elementary current amplitudes to currents.
- Return type:
ndarray
- swarmpal.toolboxes.dsecs.dsecs_algorithms.SECS_2D_CurlFree_antisym_magnetic(thetaB, phiB, thetaSECS, phiSECS, rb, rsecs, limitangle)[source]#
Calculates the mapping from antisymmetric dipolar CF SECS to magnetic field.
- Parameters:
thetaB (ndarray) – Theta coordinate of the points where magnetic field is calculated.
phiB (ndarray) – Phi coordinate of the points where magnetic field is calculated.
thetaSECS (ndarray) – Theta coordinate of the SECS poles.
phiSECS (ndarray) – Phi coordinate of the SECS poles.
rb (ndarray) – Geocentric radius of the points where the magnetic field is calculated.
rsecs (float) – Assumed radius of the spherical shell where the currents flow (km).
limitangle (ndarray) – Half-width of the uniformly distributed SECS, [radian], scalar or Nsecs-dimensional vector.
- Returns:
ndarray – 2D array that maps the SECS amplitudes to the theta component of the current density.
ndarray – 2D array that maps the SECS amplitudes to the phi component of the current density.
- swarmpal.toolboxes.dsecs.dsecs_algorithms.SECS_2D_CurlFree_antisym_vector(thetaV, phiV, thetaSECS, phiSECS, radius, limitangle)[source]#
Calculates the mapping from antisymmetric dipolar CF SECS to current density.
- Parameters:
thetaV (ndarray) – theta coordinate of points where the vector is calculated
phiV (ndarray) – phi coordinate of points where the vector is calculated
thetaSECS (ndarray) – theta coordinates od SECS poles
phiSECS (ndarray) – phi coordinates of points where the vector is calculated
radius (ndarray) – Radius of the sphere where the calculation takes place, [km], scalar
limitangle (ndarray) – Half-width of the uniformly distributed SECS, [radian], scalar or Nsecs-dimensional vector.
- Returns:
ndarray – 2D array that maps the SECS amplitudes to the theta component of the current density.
ndarray – 2D array that maps the SECS amplitudes to the phi component of the current density.
- swarmpal.toolboxes.dsecs.dsecs_algorithms.SECS_2D_DivFree_magnetic(thetaB, phiB, thetaSECS, phiSECS, rb, rsecs)[source]#
Calculates the array that maps the SECS amplitudes to magnetic field.
- Parameters:
thetaB (ndarray) – Theta coordinate of the points where magnetic field is calculated.
phiB (ndarray) – Phi coordinate of the points where magnetic field is calculated.
thetaSECS (ndarray) – Theta coordinate of the SECS poles.
phiSECS (ndarray) – Phi coordinate of the SECS poles.
rb (ndarray) – Geocentric radius of the points where the magnetic field is calculated.
rsecs (float) – Assumed radius of the spherical shell where the currents flow (km).
- Returns:
ndarray – 2D array to map the elementary current system amplitudes to radial component of magnetic field.
ndarray – 2D array to map the elementary current system amplitudes to theta component of magnetic field.
ndarray – 2D array to map the elementary current system amplitudes to phi component of magnetic field.
- swarmpal.toolboxes.dsecs.dsecs_algorithms.SECS_2D_DivFree_vector(thetaV, phiV, thetaSECS, phiSECS, radius, limitangle)[source]#
Calculates 2D DF SECS matrices for currents densities.
Function for calculating matrices matVtheta and matVphi which give the theta- and phi-components of a vector field from the scaling factors of div-free spherical elementary surrent systems (DF SECS).
- Parameters:
thetaV (ndarray) – theta coordinate of points where the vector is calculated
phiV (ndarray) – phi coordinate of points where the vector is calculated
thetaSECS (ndarray) – theta coordinates od SECS poles
phiSECS (ndarray) – phi coordinates of points where the vector is calculated
radius (ndarray) – Radius of the sphere where the calculation takes place, [km], scalar
limitangle (ndarray) – Half-width of the uniformly distributed SECS, [radian], scalar or Nsecs-dimensional vector.
- Returns:
ndarray – 2D array to map the elementary current system amplitudes to theta component of current density.
ndarray – 2D array to map the elementary current system amplitudes to phi component of current density.
- class swarmpal.toolboxes.dsecs.dsecs_algorithms.dsecsdata[source]#
Bases:
object
Class for DSECS variables and fitting procedures
- class swarmpal.toolboxes.dsecs.dsecs_algorithms.dsecsgrid[source]#
Bases:
object
Class for all the grids needed in DSECS analysis
- swarmpal.toolboxes.dsecs.dsecs_algorithms.getLocalDipoleFPtrack1D(latB, rB, Ri)[source]#
Finds the local dipole footpoints for the 1D curl-free grid creation.
- Parameters:
latB (ndarrar) – Magnetic latitude of the satellite, [degree].
rB (ndarray) – Geocentric radius of the satellite, [km].
Ri (int or float) – Assumed radius of the spherical shell where the currents flow, [km].
- Returns:
track – Latitude of the local dipole footpoints of the satellite, [degree].
- Return type:
ndarray
- swarmpal.toolboxes.dsecs.dsecs_algorithms.getLocalDipoleFPtrack2D(latB, lonB, rB, Ri)[source]#
Finds the local dipole footpoints for the 2D curl-free grid creation.
- Parameters:
latB (ndarray) – Magnetic latitude of the satellite, [degree].
lonB (ndarray) – Magnetic longitude of the satellite, [degree].
rB (ndarray) – Geocentric radius of the satellite, [km].
Ri (int or float) – Assumed radius of the spherical shell where the currents flow, [km].
- Returns:
latFP, lonFP – Latidue and longitude of the local dipole footpoints of the satellite, [degree].
- Return type:
ndarray
- swarmpal.toolboxes.dsecs.dsecs_algorithms.getUnitVectors(SwA, SwC)[source]#
Calculates the magnetic unit vectors.
- Parameters:
SwA (xarray) – Swarm A and C datasets.
SwC (xarray) – Swarm A and C datasets.
- Returns:
SwA, SwC – Swarm A and C datasets including magnetic unit vectors ‘ggUvT’, ‘ggUvP’ and ‘UvR’.
- Return type:
xarray
- class swarmpal.toolboxes.dsecs.dsecs_algorithms.grid1D[source]#
Bases:
object
Simple class to hold a 1D lat grid
- create(lat1, lat2, dlat, extLat)[source]#
Initializes the 1D grid from Swarm data.
- Parameters:
lat1 (ndarray) – Swarm A and C latitudes ,[degree].
lat2 (ndarray) – Swarm A and C latitudes ,[degree].
dlat (int or float) – 1D grid spacing in latitudinal direction, [degree].
extLat (int) – Number of points to extend the 1D grid outside data area in latitudinal direction.
- class swarmpal.toolboxes.dsecs.dsecs_algorithms.grid2D(origin='geo')[source]#
Bases:
object
Class for 2D DSECS grids
- create(lat1, lon1, lat2, lon2, dlat, lonRat, extLat, extLon, poleLat, poleLon)[source]#
Initializes the 2D grid from Swarm data.
- Parameters:
lat1 (ndarray) – Swarm A and C latitudes, [degree].
lat2 (ndarray) – Swarm A and C latitudes, [degree].
lon1 (ndarray) – Swarm A and C longitudes, [degree].
lon2 (ndarray) – Swarm A and C longitudes, [degree].
dlat (float) – 2D grid spacing in latitudinal direction, [degree].
lonRat (int or float) – Ratio between the satellite separation and longitudinal spacing of the 2D grid.
extLat (int) – Number of points to extend the 2D grid outside data area in latitudinal and longitudinal directions.
extLon (int) – Number of points to extend the 2D grid outside data area in latitudinal and longitudinal directions.
- swarmpal.toolboxes.dsecs.dsecs_algorithms.mag_transform_dsecs(SwA, SwC, pole_lat, pole_lon)[source]#
Rotates the data to a magnetic coordinate systems.
- Parameters:
SwA (xarray) – Swarm A and C datasets.
SwC (xarray) – Swarm A and C datasets.
pole_lat (float) – Latitude and longitude of the new pole in the old coordinates, [degree].
pole_lon (float) – Latitude and longitude of the new pole in the old coordinates, [degree].
- Returns:
SwA, SwC – Swarm A and C datasets including data rotated to the magnetic coordinate system.
- Return type:
xarray
- swarmpal.toolboxes.dsecs.dsecs_algorithms.secs_2d_curlFree_antisym_lineintegral(thetaB, phiB, thetaSECS, phiSECS, rb, rsecs, smallr)[source]#
Line integral for DSECS 2d curl free.
- Parameters:
thetaB (ndarray) – Theta coordinate of the points where magnetic field is calculated.
phiB (ndarray) – Phi coordinate of the points where magnetic field is calculated.
thetaSECS (float) – Theta coordinate of the SECS pole in radian.
phiSECS (float) – Phi coordinate of the SECS pole in radian.
rb (ndarray) – Geocentric radius of the points where the magnetic field is calculated.
rsecs (float) – Assumed radius of the spherical shell where the currents flow (km).
smallr (ndarray) – Limit of “small” radius [km]
- Returns:
bx, by, bz – Cartesian components of the magnetic field at the magnetometer locations [nT]
- Return type:
ndarray
- swarmpal.toolboxes.dsecs.dsecs_algorithms.trim_data(SwA, SwC)[source]#
Finds a period with suitable spaceraft velocity for analysis.
- Parameters:
SwA (xarray) – Swarm A and C datasets.
SwC (xarray) – Swarm A and C datasets.
- Returns:
SwA, SwC – Swarm A and C datasets trimmed to the suitable period.
- Return type:
xarray
swarmpal.toolboxes.dsecs.aux_tools#
- swarmpal.toolboxes.dsecs.aux_tools.get_eq(ds, lat_filter_max=60, ovals=None)[source]#
Splits data into a list of pieces suitable for DSECS analysis based on latitude.
- Parameters:
ds (xarray) – Input dataset.
lat_filter_max (int, optional) – Maximum geographic latitude (absolute value). All data above this limit is ignored. Default: 60
ovals (list) – list of pre-calculated numpy slices with which to split the data.
- Returns:
out – List of data segments split for DSECS analysis if they have been precalculated.
- Return type:
list of xarray, or None, None by default
- swarmpal.toolboxes.dsecs.aux_tools.sph2sph(latP, lonP, latOld, lonOld, VtOld, VpOld)[source]#
- Change coordinates and (horizontal part of) vectors from one spherical
coordinate system to another.
- Parameters:
latP (float) – Latitude and longitude of the new pole in the old coordinates, [deg]
lonP (float) – Latitude and longitude of the new pole in the old coordinates, [deg]
latOld (array) – Latitudes and longitudes of the data points to be transformed in the old coordinate system, [deg]
lonOld (array) – Latitudes and longitudes of the data points to be transformed in the old coordinate system, [deg]
VtOld (array) – Theta- and phi-components of a vector field at locations (latOld,lonOld). If you want to transform only coordinates make these empty arrays (np.array([])).
VpOld (array) – Theta- and phi-components of a vector field at locations (latOld,lonOld). If you want to transform only coordinates make these empty arrays (np.array([])).
- Returns:
latNew, lonNew (array) – Coordinates of the points (latOld,lonOld) in the new system, [deg], same size as input parameters. Note that the old pole has longitude 0 and longitude is between 0 and 360.
VtNew, VpNew (array) – Theta- and phi-components of the vector field in the new system. If no input vector field was given, these are empty vectors [].
- swarmpal.toolboxes.dsecs.aux_tools.sub_FindLongestNonZero(x)[source]#
- Find the longest sequence of consecutive non-zero elements in an array.
Return the sequence and indices as y=x[ind].
- Parameters:
x (array) – Input array
- Returns:
y (array) – Longest sequence of consecutive non-zero elements
ind (array) – Indices of the sequence
- swarmpal.toolboxes.dsecs.aux_tools.sub_LonInterp(lat, lon, intLat, method, ekstra=nan)[source]#
Robust interpolation of longitude by interpolating the sine and cosine of the angle. This way a 360 degree jump does not matter.
- Parameters:
lat (ndarray) – Original coordinates, [degree].
lon (ndarray) – Original coordinates, [degree].
intLat (ndarray) – Latitudes where original longitudes are interpolated to, [degree].
method (str) – Interpolation method. Same as for scipy.interpolate.interp1d.
- Returns:
intLon – Interpolated longitudes at latitudes intLat, [degree]. -180 <= intlon <= 180.
- Return type:
ndarray
- swarmpal.toolboxes.dsecs.aux_tools.sub_Swarm_grids(lat1, lon1, lat2, lon2, Dlat2D, LonRatio, ExtLat2D, ExtLon2D)[source]#
Creates 2D grids around the two Swarm satellite paths.
- Parameters:
lat1 (ndarray) – Geographic latitudes of the satellites’ paths, [degree].
lat2 (ndarray) – Geographic latitudes of the satellites’ paths, [degree].
lon1 (ndarray) – Geographic longitudes of the satellites’ paths, [degree].
lon2 (ndarray) – Geographic longitudes of the satellites’ paths, [degree].
Dlat2D (float) – 2D grid spacing in latitudinal direction, [degree].
LonRatio (float) – Ratio between the satellite separation and longitudinal spacing of the 2D grid.
ExtLat2D (int) – Number of points to extend the 2D grid outside the data area in latitudinal and longitudinal directions.
ExtLon2D (int) – Number of points to extend the 2D grid outside the data area in latitudinal and longitudinal directions.
- Returns:
lat2D, lon2D (ndarray) – Geographic latitudes and longitudes of the 2D grid, [degree].
angle2D (ndarray) – Half-angle of such a spherical cap that has the same area as the 2D grid cell, [radian].
mat2DsecondLat (ndarray) – Second gradient matrix in latitude.
- swarmpal.toolboxes.dsecs.aux_tools.sub_Swarm_grids_1D(lat1, lat2, Dlat1D, ExtLat1D)[source]#
Creates 1D grid around the two Swarm satellite paths.
- Parameters:
lat1 (ndarray) – Geographic latitudes of the satellites’ paths, [degree].
lat2 (ndarray) – Geographic latitudes of the satellites’ paths, [degree].
Dlat1D (int or float) – 1D grid spacing in latitudinal direction, [degree].
ExtLat1D (int) – Number of points to extend the 1D grid outside the data area in latitudinal direction.
- Returns:
lat1D (ndarray) – Geographic latitudes of the 1D grid, [degree].
mat1Dsecond (ndarray) – Second gradient matrix.
- swarmpal.toolboxes.dsecs.aux_tools.sub_inversion(secsMat, regMat, epsSVD, alpha, magVec)[source]#
- Solve matrix equation sysMat*resultVec = dataVec, using truncated
singular value decomposition.
- Parameters:
secsMat (ndarray) – Matrix giving the magnetic field from SECS amplitudes
regMat (ndarray) – Regularization matrix for the SECS amplitudes
epsSVD (float) – Pre-defined truncation parameter so that all singular values smaller than epsSVD*(largest singular value) will be ignored
alpha (float) – Parameter scaling the amount of regularization
magVec (array) – Magnetic measurements to be fitted
- Returns:
result_Vec – Vector of SECS amplitudes
- Return type:
array
- swarmpal.toolboxes.dsecs.aux_tools.sub_points_along_fieldline(thetaSECS, Rsecs, L, minD)[source]#
- Calculate points to be used in the Biot-Savart integral along field line
in function SECS_2D_CurlFree_AntiSym_magnetic_lineintegral
- Parameters:
thetaSECS (float) – Co-latitude of the CF SECS, [radian]
Rsecs (float) – Radius of the sphere where CF SECS is located, [km]
L (float) – L-value of the field line starting from the CF SECS pole, [km NOTE: this is in kilometers, so really L*Rsecs;
minD (float) – minimum horizontal distance between the CF SECS pole and footpoints of the points where magnetic field is needed, [km], SCALAR
- Returns:
co-latitudes of the integration points (= end points of the current elements), [radian]
- Return type:
array