TimeseriesTools.closeneighboursMethod
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.

source
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 to pointprocess!.
source
TimeseriesTools.gammarenewalMethod
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!.

source
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 to Random.default_rng().
  • t0: (optional) The initial time. Defaults to 0.0.
source
TimeseriesTools.stoicMethod
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 is npi, the integral of two gaussians with equal variance at a given distance from each other.
  • σ: Width parameter of the kernel. For npi, this is the width of the unit-mass Gaussian kernels. Default is 0.025.
  • Δt: Time window for considering spikes as close. Default is σ * 10.
source
TimeseriesTools.sttcMethod
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)
source
TimeseriesTools._energyspectrumMethod
_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.

source
TimeseriesTools._energyspectrumMethod
_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
source
TimeseriesTools._powerspectrumMethod
_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.

source
TimeseriesTools.colorednoiseMethod
colorednoise(ts::AbstractRange; α=2.0)

Generate a colored-noise time series with a specified power-law exponent α on the given times ts.

Arguments

  • ts: An AbstractRange 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
source
TimeseriesTools.energyspectrumMethod
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.

source
TimeseriesTools.powerspectrumMethod
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.

source
GeometryBasics.decomposeMethod
decompose(x::Union{<:AbstractTimeSeries, <:AbstractSpectrum})

Convert a time series or spectrum to a tuple of the dimensions and the data (as Arrays).

source
TimeseriesTools.UnitfulTimeSeriesType
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
source
TimeseriesTools.TimeSeriesMethod
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}
source
TimeseriesTools.convertconstMethod
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)
source
TimeseriesTools.dimunitMethod
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"
source
TimeseriesTools.frequnitMethod
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"
source
TimeseriesTools.timeunitMethod
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"
source
Unitful.unitMethod
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"
source