API Reference

Data Model

SpaceDataModel.AbstractDataVariableType

A variable v of a type derived from AbstractDataVariable should at least implement:

  • Base.parent(v): the parent array of the variable

Optional:

  • times(v): the timestamps of the variable
  • units(v): the units of the variable
  • meta(v): the metadata of the variable
  • name(v): the name of the variable
  • dim(v, i): the i-th dimension of the variable
  • dim(v, name): the dimension named name of the variable
source
SpaceDataModel.InstrumentType
Instrument <: AbstractInstrument

Fields

  • name: The name of the instrument
  • metadata: Additional metadata
  • datasets: Collection of datasets
source
SpaceDataModel.LDataSetType
LDataSet <: AbstractDataSet

A template for generating datasets with parameterized naming patterns.

Fields

  • format: Format string pattern for the dataset name
  • data: Dictionary of variable patterns
  • metadata: Additional metadata

Examples

using SPEDAS.MMS

# Access FPI dataset specification
lds = mms.datasets.fpi_moms

# Create a concrete dataset with specific parameters
ds = DataSet(lds; probe=1, data_rate="fast", data_type="des")

The format string and variable patterns use placeholders like {probe}, {data_rate}, which are replaced with actual values when creating a concrete DataSet.

source
SpaceDataModel.NoMetadataType
NoMetadata

Indicates an object has no metadata. But unlike using nothing, get, keys and haskey will still work on it, get always returning the fallback argument. keys returns () while haskey always returns false.

source
SpaceDataModel.ProjectType
Project <: AbstractProject

A representation of a project or mission containing instruments and datasets.

Fields

  • name: The name of the project
  • metadata: Additional metadata
  • instruments: Collection of instruments
  • datasets: Collection of datasets
source
SpaceDataModel._getfieldFunction
_getfield(v, name, default)

Return the field from a composite v for the given name, or the given default if no field is present.

See also: getfield.

source
SpaceDataModel.dimFunction
dim(x, i)
dim(x, name)

Get the i-th dimension of object x, or the dimension associated with the given name. The default implementation returns axes(x, i).

A dimension may be time-varying or dependent on other dimensions; in such cases, the effective size of the corresponding array dimension ndims can be greater than 1.

source
SpaceDataModel.getcsysMethod
getcsys(x)

Get the coordinate system of x.

If x is a instance of AbstractCoordinateSystem, return x itself. If x is a type of AbstractCoordinateSystem, return an instance of the coordinate system, i.e. x().

This is a generic function, packages should extend it for their own types.

source
SpaceDataModel.getmetaFunction
getmeta(x, key, default=nothing)

Get metadata value associated with key for object x, or default if key is not present.

source
SpaceDataModel.setmetaFunction
setmeta(x, key => value, ...; symbolkey => value2, ...)
setmeta(x, dict::AbstractDict)

Update metadata for object x for key key to have value value and return x.

source
SpaceDataModel.setmeta!Function
setmeta!(x, key => value, ...; symbolkey => value2, ...)
setmeta!(x, dict::AbstractDict)

Update metadata for object x in-place and return x. The metadata container must be mutable.

The arguments could be multiple key-value pairs or a dictionary of metadata; keyword arguments are also accepted.

Examples

setmeta!(x, :units => "m/s", :source => "sensor")
setmeta!(x, Dict(:units => "m/s", :quality => "good"))
setmeta!(x; units="m/s", calibrated=true)

Throws an error if the metadata is not mutable. Use setmeta for immutable metadata.

source

Coordinate Transformations

See GeoCotrans.jl and IRBEM.jl for more details.

SPEDAS.cotransFunction
cotrans(A, in, out; backend=:auto)
cotrans(A, out; in=get_coord(A))

Transform the data to the out coordinate system from the in coordinate system.

If backend is set to :auto (default), this function automatically chooses between Julia's GeoCotrans (if available) and Fortran's IRBEM implementation. Otherwise, it uses the specified backend.

References:

source

SPEDAS

SPEDAS.amapMethod
amap(f, a, b)

Apply a function f to the intersection of a and b.

https://github.com/rafaqz/DimensionalData.jl/issues/914

source
SPEDAS.cotransMethod
cotrans(A, in, out; backend=:auto)
cotrans(A, out; in=get_coord(A))

Transform the data to the out coordinate system from the in coordinate system.

If backend is set to :auto (default), this function automatically chooses between Julia's GeoCotrans (if available) and Fortran's IRBEM implementation. Otherwise, it uses the specified backend.

References:

source
SPEDAS.current_densityMethod
current_density(B, V)

Calculate the current density time series from the magnetic field (B) and plasma velocity (V) time series.

Assume 1-D structure along the z-direction. Remember to transform the coordinates of B and V first (e.g. using mva

source
SPEDAS.fac_matMethod
fac_mat(vec::AbstractVector; xref=[1.0, 0.0, 0.0])

Generates a field-aligned coordinate (FAC) transformation matrix for a vector.

Arguments

  • vec: A 3-element vector representing the magnetic field
source
SPEDAS.fill_gapsMethod
fill_gaps(times, data; resolution, margin)

Given a sorted vector of time stamps times and corresponding data values, this function inserts missing time stamps with a value of NaN if the gap between consecutive time stamps is larger than resolution + margin.

  • If the gap is only slightly larger (within margin of the resolution), no gap is inserted.
  • The function supports numeric times or DateTime (with appropriate resolution types).

Arguments

  • times: Sorted vector of time stamps.
  • resolution: The expected time difference between consecutive time stamps.
  • margin: Allowed deviation from resolution before inserting missing time stamps.

Returns

A tuple (full_times, full_values) where:

  • full_times is a vector containing all time stamps (original and inserted).
  • full_values is a vector of data values with NaN for inserted gaps.

References

  • https://pyspedas.readthedocs.io/en/latest/modules/pytplot/tplotmath/degap.html
source
SPEDAS.interpolate_nansMethod
interpolate_nans(da; interp=LinearInterpolation)

Interpolate only the NaN values in da along the specified dimension dims. Non-NaN values are preserved exactly as they are.

The default interpolation method interp is LinearInterpolation.

source
SPEDAS.jparallelMethod
jparallel(𝐁, curl𝐁)

Calculate the parallel component of current density with respect to magnetic field, given 𝐁 and Curl of magnetic field vector curl𝐁.

source
SPEDAS.nt2dsMethod
nt2ds(nt_arr, dim; fields=propertynames(first(nt_arr)))

Convert a NamedTuple array to a DimStack of DimArrays.

source
SPEDAS.pspectrumMethod
pspectrum(x::AbstractDimArray, spec::Spectrogram)
pspectrum(x::AbstractDimArray; nfft=256, noverlap=128, window=hamming)

Compute the power spectrum (time-frequency representation) of a time series using the short-time Fourier transform.

Returns a DimArray with frequency and original time dimensions.

See also: DSP.Spectrogram, DSP.stft

Reference

source
SPEDAS.resampleMethod
resample(arr, n; dim=1, verbose=false)

Resample an array along the dimension dim to n points. If the original length is less than or equal to n, the original array is returned unchanged.

source
SPEDAS.rotateMethod
rotate(ts::AbstractMatrix, mat::AbstractMatrix)

Coordinate-aware transformation of vector/matrix by rotation matrix(s) mat(s). Assume ts is a matrix of shape (n, 3).

source
SPEDAS.set_coordMethod

Set the coordinate system.

Updates legend names and axis labels if they include the coordinate system. Also updates the dimension name if it contains the coordinate system.

Reference:

  • https://pyspedas.readthedocs.io/en/latest/modules/pytplot/dataattgetterssetters.html#set_coords
source
SPEDAS.tinterpMethod
tinterp(A, t; interp=LinearInterpolation)

Interpolate time series A at time point(s) t using interp (default: LinearInterpolation) method. Returns interpolated value for single time point or DimArray for multiple time points.

See DataInterpolations.jl for available interpolation methods.

Examples

# Interpolate at a single time point
tinterp(time_series, DateTime("2023-01-01T12:00:00"))

# Interpolate at multiple time points using cubic spline interpolation
new_times = DateTime("2023-01-01"):Hour(1):DateTime("2023-01-02")
tinterp(time_series, new_times; interp = CubicSpline)
source
SPEDAS.tinterp_nansMethod
tinterp_nans(da::AbstractDimArray; query=timeDimType, kwargs...)

Interpolate only the NaN values in da along the specified dimensions query. Non-NaN values are preserved exactly as they are.

See also interpolate_nans

source
SPEDAS.tlingradestMethod
tlingradest(fields, positions)

Interpolate and Compute spatial derivatives such as grad, div, curl and curvature using reciprocal vector technique.

source
SPEDAS.tresampleMethod
tresample(da::DimArray, n; dim = nothing, query=nothing)

Resample a DimArray specifically along its dimension dim or query to n points. Throws an error if no dimension of type dimtype is found in the array.

source
SPEDAS.tstackMethod
tstack(vectors::AbstractVector{<:AbstractVector{T}})

Stack a time series of vectors into a matrix.

By default, each row in the output matrix represents a time point from the input vector of vectors.

source
SPEDAS.tsyncMethod
tsync(A, Bs...)

Synchronize multiple time series to have the same time points.

This function aligns the time series Bs... to match the time points of A by:

  1. Finding the common time range between all input time series
  2. Extracting the subset of A within this common range
  3. Interpolating each series in Bs... to match the time points of the subset of A

Returns a tuple containing the synchronized time series, with the first element being the subset of A and subsequent elements being the interpolated versions of Bs....

Examples

A_sync, B_sync, C_sync = tsync(A, B, C)

See also: tinterp, common_timerange

source
SPEDAS.ω2fMethod

Convert angular frequency to frequency

Reference: https://www.wikiwand.com/en/articles/Angular_frequency

source
SPEDAS.@load_project_configMacro
@load_project_config(file)

Load configuration from a file and export all key-value pairs as constants. The macro evaluates in the calling module's context.

source

Timeseries Utilities

TimeseriesUtilities.TimeseriesUtilitiesModule
TimeseriesUtilities

A collection of utilities to simplify common time series analysis.

From data cleaning to arithmetic operations (e.g. linear algebra) to common time series operations (e.g. resampling, filtering).

Data Cleaning

Query

(Windowed) Statistics

Algebra

Time-Domain Operations

Time-Frequency Domain Operations

source
TimeseriesUtilities.DiffQType
DiffQ(v, t; dim=1)

Difference quotient of v with respect to t.

To avoid undefined behavior for division by Date/DateTime, we convert the time difference to a Unitful.Quantity if eltype(v) is not a Unitful.Quantity.

source
TimeseriesUtilities.dimnumFunction
dimnum(x, query)

Get the number(s) of Dimension(s) as ordered in the dimensions of an object.

Extend the function for custom type x. By default, we fall back to DimensionalData.dimnum.

source
TimeseriesUtilities.dropnaMethod
dropna(A; dim=nothing)
dropna(A::AbstractDimArray; dim=nothing, query=nothing)

Remove slices containing NaN values along along the dim dimension.

source
TimeseriesUtilities.find_outliersMethod
find_outliers(A, [method, window]; dim = 1, kw...)

Find outliers in data A along the specified dim dimension.

Returns a Boolean array whose elements are true when an outlier is detected in the corresponding element of A.

The default method is :median (other option is :mean), which uses the median absolute deviation (MAD) to detect outliers. When the length of A is greater than 256, it uses a moving window of size 16.

See also: find_outliers_median, find_outliers_mean, isoutlier - MATLAB

source
TimeseriesUtilities.replace_outliers!Method
replace_outliers!(A, method, [find_method, window]; kwargs...)
replace_outliers!(A, method, outliers; kwargs...)

Replaces outliers in A with values determined by the specified method.

Outliers can be detected using find_outliers with optional find_method and window parameters or specified directly as a Boolean array outliers.

method can be one of the following:

  • :linear: Linear interpolation of neighboring, nonoutlier values
  • :previous: Previous nonoutlier value
  • :next: Next nonoutlier value
  • :nearest: Nearest nonoutlier value

See also: filloutliers - MATLAB

source
TimeseriesUtilities.smoothMethod
smooth(da::AbstractDimArray, window; dim=Ti, suffix="_smoothed", kwargs...)

Smooths a time series by computing a moving average over a sliding window.

The size of the sliding window can be either:

  • Quantity: A time duration that will be converted to number of samples based on data resolution
  • Integer: Number of samples directly

Arguments

  • dims=Ti: Dimension along which to perform smoothing (default: time dimension)
  • suffix="_smoothed": Suffix to append to the variable name in output
  • kwargs...: Additional arguments passed to RollingWindowArrays.rolling
source
TimeseriesUtilities.tclipMethod
tclip(A, t0, t1; query=nothing, sort=false)

Clip a Dimension or DimArray A to a time range [t0, t1].

For unordered dimensions, the dimension should be sorted before clipping (see tsort).

source
TimeseriesUtilities.tclipsMethod
tclips(xs...; trange=common_timerange(xs...))

Clip multiple arrays to a common time range trange.

If trange is not provided, automatically finds the common time range across all input arrays.

source
TimeseriesUtilities.tcrossMethod
tcross(x, y; dim = nothing, query=nothing)

Compute the cross product of two (arrays of) vectors along the specified dimension dim or query.

References:

  • https://docs.xarray.dev/en/stable/generated/xarray.cross.html
source
TimeseriesUtilities.tfilterFunction
tfilter(da, Wn1, Wn2=samplingrate(da) / 2; designmethod=nothing)

By default, the max frequency corresponding to the Nyquist frequency is used.

References

  • https://docs.juliadsp.org/stable/filters/
  • https://www.mathworks.com/help/signal/ref/filtfilt.html
  • https://docs.scipy.org/doc/scipy/reference/generated/scipy.signal.filtfilt.html

Issues

  • DSP.jl and Unitful.jl: https://github.com/JuliaDSP/DSP.jl/issues/431
source
TimeseriesUtilities.time_gridMethod
time_grid(x, dt)

Create a time grid from the minimum to maximum time in x with the step size dt.

Examples

# Create hourly time grid
time_grid(x, Hour(1))
time_grid(x, 1u"hr")

# Create 1-s intervals
time_grid(x, Second(1))
time_grid(x, 1u"second")
time_grid(x, 1u"Hz")
source
TimeseriesUtilities.timerangeMethod
timerange(times)
timerange(x1, xs...)

Get the time range (minimum and maximum) of time series data.

For a single argument, returns a tuple (tmin, tmax) containing the minimum and maximum times. For multiple arguments, returns the common time range (intersection) across all arrays - equivalent to common_timerange(x1, xs...).

Examples

# Single time series
times = [1, 2, 3, 4, 5]
timerange(times)  # (1, 5)

# Multiple time series - find common range
x1_times = [1, 2, 3, 4]
x2_times = [2, 3, 4, 5]
timerange(x1_times, x2_times)  # (2, 4)

See also: common_timerange, tminimum, tmaximum

source
TimeseriesUtilities.tmask!Method
tmask!(da, t0, t1)
tmask!(da, it::Interval)
tmask!(da, its)

Mask all data values within the specified time range(s) (t0, t1) / it / its with NaN.

source
TimeseriesUtilities.tmeanMethod
tmean(x, [dt]; dim=nothing, query=nothing)

Calculate the arithmetic mean of x along the dim dimension, optionally grouped by dt.

It returns a value if x is a vector along the dim dimension, otherwise returns a DimArray with the specified dimension dropped.

If dim is not specified, it defaults to the query dimension (dimension of type TimeDim by default).

source
TimeseriesUtilities.tmedianMethod
tmedian(x, [dt]; dim=nothing, query=nothing)

Calculate the median of x along the dim dimension, optionally grouped by dt.

It returns a value if x is a vector along the dim dimension, otherwise returns a DimArray with the specified dimension dropped.

If dim is not specified, it defaults to the query dimension (dimension of type TimeDim by default).

source
TimeseriesUtilities.tnorm_combineMethod
tnorm_combine(x; dim=nothing, name=:magnitude)

Calculate the norm of each slice along query dimension and combine it with the original components.

source
TimeseriesUtilities.toprojMethod
toproj(A, B; dim=nothing, query=nothing)

Compute vector rejection (orthogonal projection) of array A from B along specified dimension dim or query.

source
TimeseriesUtilities.tprojMethod
tproj(A, B; dim=nothing, query=nothing)

Compute vector projection of A onto B along specified dimension dim or query.

source
TimeseriesUtilities.tselectMethod
tselect(A, t, [δt]; query=nothing)

Select the value of A closest to time t within the time range [t-δt, t+δt].

Similar to DimensionalData.Dimensions.Lookups.At but choose the closest value and return missing if the time range is empty.

source
TimeseriesUtilities.tsemMethod
tsem(x, [dt]; dim=nothing, query=nothing)

Calculate the standard error of the mean of x along the dim dimension, optionally grouped by dt.

It returns a value if x is a vector along the dim dimension, otherwise returns a DimArray with the specified dimension dropped.

If dim is not specified, it defaults to the query dimension (dimension of type TimeDim by default).

source
TimeseriesUtilities.tsprojMethod
tsproj(A, B; dim=nothing, query=nothing)

Compute scalar projection of A onto B along specified dimension dim or query.

source
TimeseriesUtilities.tstdMethod
tstd(x, [dt]; dim=nothing, query=nothing)

Calculate the standard deviation of x along the dim dimension, optionally grouped by dt.

It returns a value if x is a vector along the dim dimension, otherwise returns a DimArray with the specified dimension dropped.

If dim is not specified, it defaults to the query dimension (dimension of type TimeDim by default).

source
TimeseriesUtilities.tsubtractFunction
tsubtract(x, f=nanmedian; dims=timedim(x))

Subtract a statistic (default function f: nanmedian) along dimensions (default: time dimension) from x.

source
TimeseriesUtilities.tsumMethod
tsum(x, [dt]; dim=nothing, query=nothing)

Calculate the sum of x along the dim dimension, optionally grouped by dt.

It returns a value if x is a vector along the dim dimension, otherwise returns a DimArray with the specified dimension dropped.

If dim is not specified, it defaults to the query dimension (dimension of type TimeDim by default).

source
TimeseriesUtilities.tvarMethod
tvar(x, [dt]; dim=nothing, query=nothing)

Calculate the variance of x along the dim dimension, optionally grouped by dt.

It returns a value if x is a vector along the dim dimension, otherwise returns a DimArray with the specified dimension dropped.

If dim is not specified, it defaults to the query dimension (dimension of type TimeDim by default).

source
TimeseriesUtilities.unwrapMethod
unwrap(x)

Return the innermost object of the wrapped object x with similar behavior as x (e.g. same size, same type, etc.)

source
TimeseriesUtilities.window_bf_sizesMethod
window_bf_sizes(window)

Converts a window specification to backward and forward window sizes.

When window is a positive integer scalar, the window is centered about the current element and contains window-1 neighboring elements. If window is even, then the window is centered about the current and previous elements.

source