Spatiotemporal Timeseries Prediction
An application and extension of local modeling to spatiotemporal timeseries.
Examples
Several example scripts can be found in TimeseriesPrediction/examples
. These examples are run in the examples page.
Spatio-Temporal Embeddings
some info here.
#
DelayEmbeddings.reconstruct
— Function.
reconstruct(s::AbstractArray{<:AbstractArray{T,Φ}}, em)
Reconstruct the spatial timeseries s
represented by a Vector
of AbstractArray
states using the embedding struct em
of type AbstractSpatialEmbedding
.
Returns the reconstruction in the form of a Dataset
(from DelayEmbeddings
) where each row is a reconstructed state and they are ordered first through linear indexing into each state and then incrementing in time.
#
TimeseriesPrediction.AbstractSpatialEmbedding
— Type.
AbstractSpatialEmbedding <: AbstractEmbedding
Super-type of spatiotemporal embedding methods. Valid subtypes:
SpatioTemporalEmbedding
PCAEmbedding
#
TimeseriesPrediction.SpatioTemporalEmbedding
— Type.
SpatioTemporalEmbedding{Φ,BC,X} → embedding
A spatio temporal delay coordinates structure to be used as a functor. Applies to data of Φ
spatial dimensions and gives an embedding of dimensionality X
.
embedding(rvec, s, t, α)
Operates inplace on rvec
(of length X
) and reconstructs vector from spatial timeseries s
at timestep t
and cartesian index α
. Note that there are no bounds checks for t
.
It is assumed that s
is a Vector{<:AbstractArray{T,Φ}}
.
Constructors
There are some convenience constructors that return intuitive embeddings here:
The "main" constructor is
SpatioTemporalEmbedding{X}(τ, β, bc, fsize)
which allows full control over the spatio-temporal embedding.
Χ == length(τ) == length(β)
: dimensionality of resulting reconstructed space.τ::Vector{Int}
= Vector of temporal delays for each entry of the reconstructed space (sorted in ascending order).β::Vector{CartesianIndex{Φ}}
= vector of relative indices of spatial delays for each entry of the reconstructed space.bc::BC
: boundary condition.fsize::NTuple{Φ, Int}
: Size of each state in the timeseries.
For example, if you want to have spatial neighbors [0, ±1] for time delays τ = 2, 4, then
τ = [2,2,2,4,4,4] β = CartesianIndex.([1, 0, -1, 1, 0, -1])
#
TimeseriesPrediction.cubic_shell_embedding
— Function.
cubic_shell_embedding(s, γ, τ, B, k, bc=PeriodicBoundary()) → embedding
Create a SpatioTemporalEmbedding
instance that includes spatial neighbors in hypercubic shells. The embedding is to be used with data from s
.
Description
Points are participating in the embedding by forming hypercubic shells around the current point. The total shells formed are B
. The points on the shells have spatial distance k ≥ 1
(distance in indices, like a cityblock metric). k = 1
means that all points of the shell participate. The points of the hypercubic grid can be separated by k ≥ 1
points apart (i.e. dropping k-1
in-between points). In short, in each spatial dimension of the system the cartesian offset indices are -B*k : k : k*B
.
γ
is the number of temporal steps in the past to be included in the embedding, where each step in the past has additional delay time τ::Int
. γ=0
corresponds to using only the present. Notice that all embedded time frames have the same spatial structure, in contrast to light_cone_embedding
.
As an example, consider one of the γ
embedded frames (all are the same) of a system with 2 spatial dimensions (□
= current point, (included by definition in the embedding), n
= included points in the embedding coming from n
-th shell, .
= points not included in the embedding)
B = 2, k = 1 | B = 1, k = 2 | B = 2, k = 2 | | . . . . . . . . . | . . . . . . . . . | 2 . 2 . 2 . 2 . 2 . . . . . . . . . | . . . . . . . . . | . . . . . . . . . . . 2 2 2 2 2 . . | . . 1 . 1 . 1 . . | 2 . 1 . 1 . 1 . 2 . . 2 1 1 1 2 . . | . . . . . . . . . | . . . . . . . . . . . 2 1 □ 1 2 . . | . . 1 . □ . 1 . . | 2 . 1 . □ . 1 . 2 . . 2 1 1 1 2 . . | . . . . . . . . . | . . . . . . . . . . . 2 2 2 2 2 . . | . . 1 . 1 . 1 . . | 2 . 1 . 1 . 1 . 2 . . . . . . . . . | . . . . . . . . . | . . . . . . . . . . . . . . . . . . | . . . . . . . . . | 2 . 2 . 2 . 2 . 2
#
TimeseriesPrediction.light_cone_embedding
— Function.
light_cone_embedding(s, γ, τ, r0, c, bc=PeriodicBoundary()) → embedding
Create a SpatioTemporalEmbedding
instance that includes spatial and temporal neighbors of a point based on the notion of a light cone.
The embedding is to be used with data from s
.
Description
Information does not travel instantly but with some finite speed c ≥ 0.0
. This constructor creates a cone-like embedding including all points in space and time, whose value can influence a prediction based on the information speed c
. γ
is the number of temporal steps in the past to be included in the embedding, where each step in the past has additional delay time τ::Int
. γ=0
corresponds to using only the present. r0
is the initial radius at the top of the cone, i.e. the radius of influence at the present. bc
is the boundary condition.
The radius of the light cone evolves as: r_i = i*τ*c + r0
for each step i ∈ 0:γ
.
As an example, in a one-dimensional system with γ = 1, τ = 2, r0 = 1
, the embedding looks like (□
= current point (included by definition in the embedding), o
point to be predicted using temporalprediction
, x
= points included in the embedding, .
= points not included in the embedding)
time | c = 1.0 | c = 2.0 | c = 0.0 n + 1 | ..........o.......... | ..........o.......... | ..........o.......... n | .........x□x......... | .........x□x......... | .........x□x......... n - 1 | ..................... | ..................... | ..................... n - τ | .......xxxxxxx....... | .....xxxxxxxxxx...... | .........xxx.........
Besides this example, in the official documentation we show a function explain_light_cone
which produces a plot of the light cone for 2 spatial dimensions (great for understanding!).
#
TimeseriesPrediction.PCAEmbedding
— Type.
PCAEmbedding(s, em::AbstractSpatialEmbedding; kwargs...) → embedding
A spatio temporal delay coordinates structure with Principal Component Analysis as a means of dimension reduction, embedding
can be used as a functor:
embedding(rvec, s, t, α)
which operates inplace on rvec
and reconstructs vector from spatial time series s
at timestep t
and cartesian index α
.
To instantiate this embedding
, give the data to be reconstructed s
as well as an instance of SpatioTemporalEmbedding
to PCAEmbedding
.
Keyword Arguments
pratio = 0.99
: Ratio of variances that needs to be preserved in low-dimensional PCA-reconstruction.maxoutdim = 25
: Upper limit for output dimension. May breakpratio
criterion.every_t = 1
: Speed up computation by only using every n-th point in time.every_α = 1
: Speed up computation further by only using every n-th point in space (linear indexing).
To set the output dimension to a certain value X
, pass pratio=1, maxoutdim=X
.
Here is a function that visualizes how the light_cone_embedding
works in 2 spatial dimensions:
using PyPlot, TimeseriesPrediction using LinearAlgebra: norm function explain_light_cone(;γ = 2, τ = 2, r = 1, c = 1) maxr = D*τ*c + r0 figure() xticks(-maxr:maxr) yticks(-maxr:maxr) for i in γ:-1:0 r = i*τ*c + r points = TimeseriesPrediction.indices_within_sphere(r, 2) radius = maximum(norm(Tuple(p)) for p in points) if r != 0 x = r*cos.(range(0, stop = 2π, length = 100)) y = r*sin.(range(0, stop = 2π, length = 100)) plot(x,y, c = "C$i") end x = map(xy -> xy[1], points) y = map(xy -> xy[2], points) scatter(x,y, c = "C$i", s=100, zorder = 3, label = "within \$r = r + $i \\tau c\$") end title("γ = $γ, τ = $τ, r = $r, c = $c") PyPlot.grid(zorder = -1) legend(loc = "upper right") xlabel("x") ylabel("y") PyPlot.axis("equal") tight_layout() end
Boundary conditions
#
TimeseriesPrediction.ConstantBoundary
— Type.
ConstantBoundary(b) <: AbstractBoundaryCondition
Constant boundary condition type. Enforces constant boundary conditions when passed to SpatioTemporalEmbedding
by filling missing out-of-bounds values in the reconstruction with parameter b
.
#
TimeseriesPrediction.PeriodicBoundary
— Type.
PeriodicBoundary <: AbstractBoundaryCondition
Periodic boundary condition struct. Enforces periodic boundary conditions when passed to SpatioTemporalEmbedding
in the reconstruction.
Prediction functions
#
TimeseriesPrediction.temporalprediction
— Function.
temporalprediction(U, em::AbstractSpatialEmbedding, tsteps; kwargs...)
Perform a spatio-temporal time series prediction for tsteps
iterations, using local weighted modeling [1] give a time series of the form U::AbstractVector{<:AbstractArray{T, Φ}}
.
The returned data always contains the final state of U
as starting point (total returned length is tsteps+1
). The reconstruction process is defined by em
. For available methods and interfaces see AbstractSpatialEmbedding
.
Keyword Arguments
ttype = KDTree
: Type/Constructor of tree structure. So far only tested withKDTree
.method = AverageLocalModel(ω_safe)
: Subtype ofAbstractLocalModel
.ntype = FixedMassNeighborhood(3)
: Subtype ofAbstractNeighborhood
.initial_ts = U
: Initial states for prediction (same type asU
). Must have at least as many states as the maximum delay time used. Defaults to the training setU
.progress = true
: To print progress done.
Description
This method works similarly to localmodel_tsp
, by expanding the concept of delay embedding to spatially extended systems. Instead of reconstructing complete states of the system, local states are used. See AbstractSpatialEmbedding
for details on the embedding. Predictions are then performed frame by frame and point py point. Once all values for a new frame are found, the frame is added to the end of the timeseries and used to generate new prediction queries for the next time step.
Performance Notes
Be careful when choosing embedding parameters as memory usage and computation time depend strongly on the resulting embedding dimension.
References
[1] : U. Parlitz & C. Merkwirth, Phys. Rev. Lett. 84, pp 1890 (2000)
#
TimeseriesPrediction.crossprediction
— Function.
crossprediction(source_train, target_train, source_pred, em::AbstractSpatialEmbedding; kwargs...)
Perform a spatio-temporal timeseries cross-prediction for target
from source
, using local weighted modeling [1]. This can be used for example when there are coupled spatial fields and one is used to predict the other. It is assumed that source_train
, target_train
, source_pred
are all of the same type, AbstractVector{<:AbstractArray{T, Φ}}
.
The spatio temporal delay embedding process is defined by em
. See AbstractSpatialEmbedding
for available methods and interfaces.
Keyword Arguments
ttype = KDTree
: Type/Constructor of tree structure. So far only tested withKDTree
.method = AverageLocalModel(ω_safe)
: Subtype ofAbstractLocalModel
.ntype = FixedMassNeighborhood(3)
: Subtype ofAbstractNeighborhood
.progress = true
: To print progress done.
Description
The reconstructed state of source_train[t][i,j,...]
is associated with the output value target_train[t][i,j,...]
. This establishes a "connection" between target
and source
. Taking a reconstructed state of source_pred
as query point, the function finds its neighbors in the reconstructed space of source_train
using neighborhood ntype
. Then, the neighbor indices are used to make a prediction for the corresponding value of the target
, using the established "connection" between fields.
Additional Interfaces
To save computation time in the case of repeated predictions with the same training set and embedding parameters we provide an additional interface that allows the user to provide an existing reconstruction and tree structure.
R = reconstruct(train_in, em) tree = ttype(R) params = PredictionParams(em, method, ntype, ttype) sol = crossprediction(params, train_out, pred_in, R, tree; progress=true)
where params
is an internal container with all relevant parameters.
Performance Notes
Be careful when choosing embedding parameters as memory usage and computation time depend strongly on the resulting embedding dimension.
References
[1] : U. Parlitz & C. Merkwirth, Phys. Rev. Lett. 84, pp 1890 (2000)