CyclingSignatures.curve_cycleMethod
curve_cycle(i,j;F=DEFAULT_FIELD)

Returns the curve cycle corresponding to i and j where i<j.

Returns

  • edges: vector of edges [i,i+1], [i+1,i+2], ..., [j-1,j], [i,j].
  • coeffs: coefficients [1,...,1,-1].
source
CyclingSignatures.dm_components_explicitMethod
dm_components_explicit(points, metric, flt_threshold)

Computes a filtration-minimal vertex of each connected component of the distance complex up to the filtration threshold flt_threshold. The distance complex is specified by points and metric.

This is implmented using a two pass annotation algorithm for images, the image being the upper right part of the distance matrix. In this implementation, the distance matrix is computed explicitly.

source
CyclingSignatures.dm_components_first_pass_explicitMethod
dm_components_first_pass_explicit(points, metric, threshold)

Given a matrix with trajectory points, a metric and a filtrationThreshold, computes the connected components of the graph defined by the matrix pairwise(metric,points) .<= threshold where

  • vertices are matrix entries,
  • edges are nearest 4.

Returns a disjoint set with labels, the label of each entry and the distance matrix. Note: entries with different labels may be in the same component (iff their labels are in the same set in cc)

source
CyclingSignatures.dm_components_first_pass_implicitMethod
dm_components_first_pass_explicit(points, metric, threshold)

Given a matrix with trajectory points, a metric and a filtration threshold, computes the connected components of the graph defined by the matrix pairwise(metric,points) .<= threshold where

  • vertices are matrix entries,
  • edges are nearest 4.

Returns a disjoint set with label and the label of each entry. Note: entries with different labels may be in the same component (iff their labels are in the same set in cc)

source
CyclingSignatures.dm_components_implicitMethod
filtration_min_vertices_implicit_dm(points, metric, flt_threshold)

Computes a filtration-minimal vertex of each connected component of the distance complex specified by points and metric up to the filtration threshold flt_threshold. This is implmented using a two pass annotation algorithm for images, the image being the upper right part of the distance matrix. In this implementation, the distance matrix is not computed explicitly.

source
CyclingSignatures.trajectory_barcodeFunction
trajectory_barcode(::Val{:DistanceMatrix}, points, metric, fltThreshold, field=DEFAULT_FIELD)

Computes a trajectory barcode using the distance matrix method, for arguments see trajectory_barcode. This assumes (and formall makes use of) a curve hypothesis.

The distance complex of a the points in points with respect to the metric is the filtered complex where

  • vertices are pairs (i, j) with i < j, filtered by metric(points[:, i], points[:, j]) and
  • edges are nearest 4, with filtration value of an edge being the max of its adjacent vertices.

Note:

In general, the returned collection of bars is not a persistence diagram for the point cloud.

source
CyclingSignatures.vertex_to_essential_barMethod
vertex_to_essential_bar(node, points, metric; field=DEFAULT_FIELD)

Given a node (i,j) in the distance matrix, returns an essential bar corresponding to the induced cycle in the Vietoris–Rips complex.

source
CyclingSignatures.colspace_normal_formMethod
colspace_normal_form!(M)

Computes and returns the (nonzero part of the) reduced column echelon form of the matrix M.

Elementary column operations are performed to transform M into a form where:

  • Each pivot value is 1 and is the only non-zero entry in its column.
  • The leading entry of each non-zero column is below the leading entry of the previous row.
  • columns without leading entries are removed from the matrix.
source
CyclingSignatures.AbstractCubicalAcyclicCarrierType

Abstract base type for cubical acyclic carriers. Subtypes represent an implementation of a cell complex representing a union of cubes which allows the construction of chains from acyclic subsets and the annotation of cycles using a cohomology basis.

Subtypes must implement:

  • induced_one_chain(carrier, edge_boxes): Given a vector of boxes covering an edge, computes the induced one-chain.
  • annotate_chain(carrier, chain): Annotates a chain using a cohomology basis. For cycles, the result is independent of the representative.
  • betti_1(carrier): Returns the first betti number of the carrier.
source
CyclingSignatures.AbstractCubicalComparisonSpaceType

Abstract type for cubical-based comparison spaces. Inherits from AbstractComparisonSpace.

Subtypes must implement:

  • edge_boxes(comparison_space, p1, p2): Computes the boxes covering the edge between points p1 and p2.
  • carrier(comparison_space): Returns an associated cubical acyclic carrier.
  • betti_1(comparison_space): Returns the first betti number of the comparison space.
source
CyclingSignatures.AbstractSampleableTrajectoryType

Abstract base type for trajectories that can be used in cycling computations.

Subtypes must implement:

  • evaluate_interval(trajectory, a, b): Returns points that the trajectory visits between time a and b.
  • time_domain(trajectory): Returns time domain of the trajectory
source
CyclingSignatures.CubicalVRCarrierType

Fields

  • pts::Matrix{S}: Points in the carrier, where S is an integer type. The columns are assumed to be sorted.
  • cplx::CellComplex{Simplex}: The Nerve of the cubical cover, it is the $l_\infty$-distance VR-complex of the box centers.
  • h1::V: The transpose of h1 is a basis for the first cohomology of the complex.
source
CyclingSignatures.DynamicDistanceType
DynamicDistance{C} <: Metric

Calculates the distance between points in the tangent bundle according to the formula

$d_{dyn}((p,v),(q,w)) = \max(\|p - q\|, C \|v - w\|)$

where C is a constant.

Fields

  • dim::Int: The dimension of the space (i.e. half the tangent space dimension).
  • c::C: The constant used in the distance calculation.
source
CyclingSignatures.FFType
FF{M} <: Integer

FF{M} is the default field used in CyclingSignatures. It is a representation of a finite field $\mathbb{Z}_M$, integers modulo small, prime M. Supports field arithmetic and can be converted to integer with Int.

Example

julia> FF{3}(5)
2 mod 3

julia> FF{3}(5) + 1
0 mod 3

Copyright

Copyright (c) 2020 mtsch

source
CyclingSignatures.RefinedEquidistantTrajectoryType
struct RefinedEquidistantTrajectory

Models a refinement of a trajectory on an equidistant grid.

More precisely, the points [trajectory[:, i] for i in t_vec[1:end-1]] are assumed to be the sample points, and the other points are the refined points.

It holds for all i in 1:length(t_vec)-1 that the points in the i-th interval are given by the indices in t_vec[i]:t_vec[i+1]-1.

source
CyclingSignatures.cycling_signatureFunction
cycling_signature([alg::Val, ]trajectory_space::TrajectorySpace, range; r_max, field=DEFAULT_FIELD)

Computes the cycling signature of the segment specified by range inside of trajectory_space for filtration values $[0,r_{max}]$. Optionally, a field and an algorithm can be specified.

Arguments

  • alg: currently only :DistanceMatrix works
  • trajectory_space: the trajectory space
  • range: currently, this evaluates the cycling signature on the interval [first(range):last(range)]
  • r_max: if unspecified, the default in trajectory_space is used
  • field: coefficient field for homology computation

Returns

TODO

source
CyclingSignatures.cycspace_inclusion_matrixMethod
cycspace_inclusion_matrix(V, W)

For a list of matrices V and W, compute the inclusion matrix where entry (i,j) is 1 if the image of V[i] is included in the image of W[j], and 0 otherwise.

source
CyclingSignatures.dynamicIndicesMethod
function dynamicIndices(dataq::Matrix)

Returns a vector v and a matrix M s.t. the columns of M are pairwise distinct and M[:,v[t]] = dataq[:,t] for all t. In other words, M contains the unique boxes and v[t] the index of the box at time step t.

source
CyclingSignatures.fix_edgeMethod
fix_edge(carrier::CubicalVRCarrier, box1, box2)

Searches for a (short) shortest path between the boxes. Note that current implementation is not well-defined since it may be random which way a missing box is bypassed.

source
CyclingSignatures.induced_one_chainMethod
induced_one_chain(carrier::AbstractCubicalAcyclicCarrier, edge_boxes)

Implements the conversion of an acyclic cubical cover to an induced chain.

Arguments

  • carrier::AbstractCubicalAcyclicCarrier: The cubical acyclcic carrier
  • edge_boxes: Boxes that cover an one chain

Returns

A representation of an induced chain

source
CyclingSignatures.lenght_countdict_to_countmatrixMethod
lenght_countdict_to_countmatrix(ct_dict; n_relevant = nothing)

Convert a length count dictionary (i.e. key maps to a vector of counts per segment length) to a count matrix for the n_relevant most appearing keys.

source
CyclingSignatures.map_cycleMethod
map_cycle(cs::AbstractComparisonSpace, points, simplices, coeffs)

Annotates the cycle specified by points, simplices and coeffs into a comparison space using the given comparison space.

Arguments

  • cs::AbstractComparisonSpace: The comparison space
  • points: Points defining the cycle
  • simplices: Simplices in the cycle
  • coeffs: Coefficients for the cycle

Returns

Mapped cycle representation (implementation-dependent)

source
CyclingSignatures.projection_max_switchesMethod
projection_max_switches(v1, v2)

Computes a vector of times T which partition the interval [0,1] such that λ \mapsto (1-λ)v_1 + λ v_2 has the same maximal component in each interval [T[0];T[1]].

source
CyclingSignatures.quantizeMethod

quantizes each column of x to a grid of size r.

If r is an integer, boxes have size r in each dimension, if r is a vector, it has to contain a length for every dimension.

source
CyclingSignatures.resampleToConsistentMethod
function resampleToConsistent(ds, Y, r, dt)

Adds additional sample points to a time series Y from a dynamical system ds at time step dt such that the resulting time series quantized at r is consistent.

Instead of the actual trajectory points, one can also make pp(Y) dynamically consistent, where pp is a postprocessing applied to the data (for example a rescaling).

Furthermore, it is possible to specify a sbradius and a sbfct which ensures consistency in the sphere bundle. Note: this works as follows: the output of the integrator is applied to sbfct, this will be scaled to sbRadius (in lInf norm) and then checked for consistency

source
CyclingSignatures.resampleToConsistentMethod
function resampleToConsistent(ds, ic, ec, r, dt; pp=identity, verbose=false, max_depth=16)

Returns a matrix M such that [M ec] is dynamically consistent and M[:,1] = ic. Note that M[:,end] != ec.

source
CyclingSignatures.resampleToDistanceFunction
function resampleToDistance(ds, Y, r, dt; kwargs...)

Adds additional sample points to a time series Y from a dynamical system ds at time step dt such that the resulting time series satisfies the specified distance requirements

Instead of the actual trajectory points, one can also require pp(Y) to satisfy the requirement, where pp is a postprocessing applied to the data (for example a rescaling).

Furthermore, it is possible to specify a sbradius and a sbfct which ensures consistency in the sphere bundle. Note: this works as follows: the output of the integrator is applied to sbfct, this will be scaled to sbRadius (in lInf norm) and then checked for consistency

source