TimeseriesTools.closeneighbours Method
closeneighbours(x, y; Δt)
Constructs a sparse matrix of distances between neighbouring spikes in two sorted spike trains.
Arguments
x
: A sorted array representing the first spike train.y
: A sorted array representing the second spike train.Δt
: The maximum time difference allowed for two spikes to be considered neighbours.
Returns
A sparse matrix D
where D[i, j]
represents the distance between the i
-th spike in x
and the j
-th spike in y
, for pairs of spikes within Δt
of each other.
TimeseriesTools.gammarenewal! Method
gammarenewal!(spikes, α, θ; t0 = randn() * α / θ, kwargs...)
Generate a sequence of gamma-distributed renewal spikes.
Arguments:
spikes::AbstractVector
: The vector to store the generated spikes.α
: The shape parameter of the gamma distribution.θ
: The scale parameter of the gamma distribution.t0
: The initial time of the spike train. Defaults to a random value drawn from a normal distribution with mean of 0 and standard deviation equal to the mean firing rate.kwargs...
: Additional keyword arguments to be passed topointprocess!
.
TimeseriesTools.gammarenewal Method
gammarenewal(N, α, θ; t0)
Generate a spike train with inter-spike intervals drawn from a Gamma process.
Arguments
N
: Number of spikes to generate.α
: Shape parameter of the gamma distribution (equivalent to the mean ISI divided by the
Fano factor).
θ
: Scale parameter of the gamma distribution (equivalent to the Fano factor).t0
: The initial time of the spike train, prior to the first sampled spike. Defaults to a
random value drawn from a normal distribution with mean of 0 and standard deviation equal to the mean firing rate.
Returns
- A
SpikeTrain
containing the generated spike times.
See also gammarenewal!
.
TimeseriesTools.pointprocess! Method
pointprocess!(spikes, D::Distribution; rng = Random.default_rng(), t0 = 0.0)
Simulate a point process by sampling inter-spike intervals from a given distribution.
Arguments
spikes
: An array to store the generated spike times.D::Distribution
: The distribution from which to sample inter-spike intervals.rng
: (optional) The random number generator to use. Defaults toRandom.default_rng()
.t0
: (optional) The initial time. Defaults to0.0
.
TimeseriesTools.stoic Method
stoic(a, b; kpi = npi, σ = 0.025, Δt = σ * 10)
Compute the spike-train overlap-integral coefficient between two spike trains, after normalizing both convolutions to unit energy
See the unnamed metric from "Schreiber S, Fellous JM, Whitmer JH, Tiesinga PHE, Sejnowski TJ (2003). A new correlation based measure of spike timing reliability. Neurocomputing 52:925-931."
Arguments
a
: Spike train a.b
: Spike train b.kpi
: Kernel product integral, a function of the distance between two spikes. Default isnpi
, the integral of two gaussians with equal variance at a given distance from each other.σ
: Width parameter of the kernel. Fornpi
, this is the width of the unit-mass Gaussian kernels. Default is0.025
.Δt
: Time window for considering spikes as close. Default isσ * 10
.
TimeseriesTools.sttc Method
sttc(a, b; Δt = 0.025)
The spike-time tiling coefficient, a measure of correlation between spike trains [1].
Arguments
a::Vector{<:Real}
: A sorted vector of spike times.b::Vector{<:Real}
: A second sorted vector of spike times .Δt::Real=0.025
: The time window for calculating the STTC.
Returns
sttc::Real
: The STTC value.
References
[1] [Cutts & Eglen 2014](https://doi.org/10.1523%2FJNEUROSCI.2767-14.2014)
TimeseriesTools.AbstractSpectrum Type
AbstractSpectrum{T, N, B}
A type alias for an AbstractToolsArray
in which the first dimension is 𝑓
requency.
TimeseriesTools.FreqIndex Type
FreqIndex
A type alias for a tuple of dimensions, where the first dimension is of type FrequencyDim
.
TimeseriesTools.MultivariateSpectrum Type
MultivariateSpectrum{T} = AbstractSpectrum{T, 2} where T
A type alias for a multivariate spectrum.
sourceTimeseriesTools.RegularFreqIndex Type
RegularFreqIndex
A type alias for a tuple of dimensions, where the first dimension is a regularly sampled 𝑓
requency.
TimeseriesTools.RegularSpectrum Type
RegularSpectrum{T, N, B}
A type alias for a spectrum with a regularly sampled frequency index.
sourceTimeseriesTools.UnivariateSpectrum Type
UnivariateSpectrum{T} = AbstractSpectrum{T, 1} where T
A type alias for a univariate spectrum.
sourceTimeseriesTools.Spectrum Method
Spectrum(f, v, x)
Constructs a multivariate spectrum with frequencies f
, variables v
, and data x
.
TimeseriesTools.Spectrum Method
Spectrum(f, x)
Constructs a univariate spectrum with frequencies f
and data x
.
TimeseriesTools._energyspectrum Method
_energyspectrum(x::RegularTimeSeries, f_min=samplingrate(x)/min(length(x)÷4, 1000))
Computes the energy spectrum of a regularly sampled time series x
with an optional minimum frequency f_min
.
TimeseriesTools._energyspectrum Method
_energyspectrum(x::RegularTimeSeries, f_min=0)
Computes the energy spectrum of a time series using the fast Fourier transform.
If f_min > 0
, the energy spectrum is calculated for windows of the time series determined by f_min
, the minimum frequency that will be resolved in the spectrum. If f_min > 0
, the second dimension of the output will correspond to the windows. For an averaged periodogram, see energyspectrum
.
If the input time series is a UnitfulTimeSeries
, the frequency will also have units. Moreover if the elements of x
are unitful, so are the elements of the spectrum.
Examples
julia> using TimeseriesTools
julia> t = range(0.0, stop=1.0, length=1000);
julia> x = sin.(2 * π * 5 * t);
julia> ts = RegularTimeSeries(t, x);
julia> S = _energyspectrum(ts);
julia> S isa MultivariateSpectrum
TimeseriesTools._powerspectrum Method
_powerspectrum(x::AbstractTimeSeries, f_min=samplingrate(x)/min(length(x)÷4, 1000); kwargs...)
Computes the power spectrum of a time series x
in Welch periodogram windows. Note that the _powerspectrum
is simply the _energyspectrum
divided by the duration of each window. See _energyspectrum
.
TimeseriesTools.colorednoise Method
colorednoise(ts::AbstractRange; α=2.0)
Generate a colored-noise time series with a specified power-law exponent α
on the given times ts
.
Arguments
ts
: AnAbstractRange
representing the time range of the generated noise.α
: The power-law exponent of the colored noise, which will have a spectrum given by 1/f^α. Defaults to 2.0.
Returns
- A
TimeSeries
containing the generated colored noise.
Example
julia> using TimeseriesTools
julia> pink_noise = colorednoise(1:0.01:10; α=1.0)
julia> pink_noise isa RegularTimeSeries
TimeseriesTools.energyspectrum Method
energyspectrum(x::RegularTimeSeries, f_min=0; kwargs...)
Computes the average energy spectrum of a regularly sampled time series x
. f_min
determines the minimum frequency that will be resolved in the spectrum. See _energyspectrum
.
TimeseriesTools.powerspectrum Method
powerspectrum(x::AbstractTimeSeries, f_min=samplingrate(x)/min(length(x)÷4, 1000); kwargs...)
Computes the average power spectrum of a time series x
using the Welch periodogram method.
TimeseriesTools.UnitfulIndex Type
UnitfulIndex
A type alias for a union of AbstractArray
, AbstractRange
, and Tuple
types with Unitful.Time
elements.
TimeseriesTools.UnitfulSpectrum Type
UnitfulSpectrum{T,N,B}
A type representing spectra with unitful frequency units.
sourceTimeseriesTools.UnitfulTimeIndex Type
UnitfulTimeIndex
A type alias for a tuple of dimensions, where the first dimension is of type DimensionalData.Dimension{<:UnitfulIndex}
.
TimeseriesTools.UnitfulTimeSeries Type
UnitfulTimeSeries{T, N, B}
A type alias for an AbstractToolsArray
with a UnitfulTimeIndex
.
Examples
julia> using Unitful;
julia> t = (1:100)u"s";
julia> x = rand(100);
julia> uts = TimeSeries(t, x);
julia> uts isa UnitfulTimeSeries
TimeseriesTools.TimeSeries Method
TimeSeries(t, x, timeunit::Unitful.Units)
Constructs a univariate time series with time t
, data x
and time units specified by timeunit
. Note that you can add units to the elements of a time series x
with, for example, x*u"V"
.
Examples
julia> using Unitful;
julia> t = 1:100;
julia> x = rand(100);
julia> ts = TimeSeries(t, x, u"ms")*u"V";
julia> ts isa Union{UnivariateTimeSeries, RegularTimeSeries, UnitfulTimeSeries}
TimeseriesTools.convertconst Method
TimeseriesTools.convertconst(c::Number, u::Unitful.Quantity)
Converts a constant c
to have the same units as u
.
Examples
julia> using Unitful;
julia> c = 5;
julia> u = 3u"s";
julia> converted_c = TimeseriesTools.convertconst(c, u);
julia> typeof(converted_c) == typeof(u)
TimeseriesTools.dimunit Method
dimunit(x::UnitfulTimeSeries, dim)
Returns the unit associated with the specified dimension dim
of a UnitfulTimeSeries
.
Examples
julia> using Unitful;
julia> t = 1:100;
julia> x = rand(100);
julia> ts = TimeSeries(t, x, u"ms");
julia> TimeseriesTools.dimunit(ts, 𝑡) == u"ms"
TimeseriesTools.frequnit Method
frequnit(x::UnitfulSpectrum)
Returns the frequency units associated with a UnitfulSpectrum
.
Examples
julia> using Unitful;
julia> t = 1:100;
julia> x = rand(100);
julia> ts = TimeSeries(t, x, u"ms");
julia> sp = fft(ts); # assuming fft returns a UnitfulSpectrum
julia> frequnits(sp) == u"Hz"
TimeseriesTools.timeunit Method
timeunit(x::UnitfulTimeSeries)
Returns the time units associated with a [UnitfulTimeSeries
].
Examples
julia> using Unitful;
julia> t = 1:100;
julia> x = rand(100);
julia> ts = TimeSeries(t, x, u"ms");
julia> timeunit(ts) == u"ms"
Unitful.unit Method
unit(x::AbstractArray)
Returns the units associated with the elements of an UnitfulTimeSeries
or UnitfulSpectrum
.
Examples
julia> using Unitful;
julia> t = 1:100;
julia> x = rand(100);
julia> ts = TimeSeries(t, x, u"ms")*u"V";
julia> unit(ts) == u"V"