Graph Laplacian/diffusion maps-based LCS methods

The LCS approaches implemented and described in this section are strongly influenced by ideas developed in the spectral clustering/diffusion maps communities. The general references include:

In the LCS context, these ideas have been adopted in the following works:

For demonstrations on example cases, please consult the page on Working with trajectories.

Function documentation

Sparsification methods

Three commonly used sparsification methods are implemented for use with various graph Laplacian methods.

CoherentStructures.KNNType
KNN(k) <: SparsificationMethod

Defines the KNN (k-nearest neighbors) sparsification method. In this approach, first k nearest neighbors are sought. In the final graph Laplacian, only those particle pairs are included which are contained in some k-Neighborhood.

source
CoherentStructures.MutualKNNType
MutualKNN(k) <: SparsificationMethod

Defines the mutual KNN (k-nearest neighbors) sparsification method. In this approach, first k nearest neighbors are sought. In the final graph Laplacian, only those particle pairs are included which are mutually contained in each others k-Neighborhood.

source
CoherentStructures.NeighborhoodType
Neighborhood(ε) <: SparsificationMethod

Defines the ε-Neighborhood sparsification method. In the final graph Laplacian, only those particle pairs are included which have distance less than ε.

source

Other sparsification methods can be implemented by defining a SparsificationMethod type and a corresponding spdist method.

Diffusion-maps-type/adjancency-matrix-based graph Laplacian methods

Since the Euclidean heat kernel is ubiquitous in diffusion maps-based computations, we provide it for convenience.

CoherentStructures.gaussianFunction
gaussian(σ=1.0)

Returns the Euclidean heat kernel as a callable function

\[x \mapsto \exp\left(-\frac{x^2}{4\sigma}\right)\]

Example

julia> kernel = gaussian(2.0);

julia> kernel(0.)
1.0
source

To compute a sparse distance matrix (or adjacency matrix, depending on the sparsification method), use spdist.

CoherentStructures.spdistFunction
spdist(data, sp_method, metric=Euclidean()) -> SparseMatrixCSC

Return a sparse distance matrix as determined by the sparsification method sp_method and metric.

source

The main high-level functions for the above-listed methods are the following.

CoherentStructures.sparse_diff_op_familyFunction
sparse_diff_op_family(data, sp_method, kernel, op_reduce; α, metric, verbose)

Return a list of sparse diffusion/Markov matrices P.

Arguments

  • data: a list of trajectories, each a list of states of type SVector;
  • sp_method: a sparsification method;
  • kernel=gaussian(): diffusion kernel, e.g., gaussian;
  • op_reduce=P -> prod(LinearMaps.LinearMap, reverse(P)): time-reduction of diffusion operators, e.g. Statistics.mean (space-time diffusion maps), P -> row_normalize!(min.(sum(P), 1)) (network-based coherence) or the default P -> prod(LinearMaps.LinearMap, reverse(P)) (time coupled diffusion maps)

Keyword arguments

  • α=1: exponent in diffusion-map normalization;
  • metric=Euclidean(): distance function w.r.t. which the kernel is computed, however, only for point pairs where $metric(x_i, x_j)\leq \varepsilon$;
  • verbose=false: whether to print intermediate progress reports.
source
CoherentStructures.sparse_diff_opFunction
sparse_diff_op(data, sp_method, kernel; α=1, metric=Euclidean()) -> SparseMatrixCSC

Return a sparse diffusion/Markov matrix P.

Arguments

  • data: a list of trajectories, each a list of states of type SVector, or a list of states of type SVector;
  • sp_method: a sparsification method;
  • kernel: diffusion kernel, e.g., gaussian;

Keyword arguments

  • α: exponent in diffusion-map normalization;
  • metric: distance function w.r.t. which the kernel is computed, however, only for point pairs where $metric(x_i, x_j)\leq \varepsilon$.
source

Normalization functions

In the diffusion maps framework, there are two commonly used normalization steps:

  1. kernel-density estimate normalization (kde_normalize!), and
  2. row normalization (row_normalize!), to obtain a diffusion/Markov operator (w.r.t. right- and left-action, respectively).
CoherentStructures.kde_normalize!Function
kde_normalize!(A, α=1)

Normalize rows and columns of A in-place with the respective row-sum to the α-th power; i.e., return $a_{ij}:=a_{ij}/q_i^{\alpha}/q_j^{\alpha}$, where $q_k = \sum_{\ell} a_{k\ell}$. Default for α is 1.

source

Diffusion-coordinate-like functions

CoherentStructures.diffusion_coordinatesFunction
diffusion_coordinates(P, n_coords) -> (Σ, Ψ)

Compute the (time-coupled) diffusion coordinate matrix Ψ and the coordinate weight vector Σ for a diffusion operator P. The argument n_coords determines the number of diffusion coordinates to be computed.

source