ddHodge

Documentation for ddHodge.

ddHodge.alignSpacesMethod
alignSpaces(graph, subspaces) -> (alignedSpaces, sheafGrad, restrictionMaps)

Solves subspace alignment problem by minimizing Rayleigh quotient of sheaf (connection) laplacian $\Delta_s$

\[\min_{\bm x} \sum_{(u,v) \in E} \|O_u \bm x_u - O_v \bm x_v\|^2 = \bm x^\top \Delta_s \bm x\]

where $O_u, O_v$ are the restriction maps of parallel transport matrix.

source
ddHodge.basischangeMethod
basischange(A, V, W) -> B::Matrix

Change of basis: the matrix A in the basis V into the basis W

Warning

Use with caution when span(V) ≠ span(W).

source
ddHodge.basischangeMethod
basischange(A, V) -> B::Matrix

Change of basis: A in the basis V into the canonical basis, i.e., W = I

source
ddHodge.ddHodgeWorkflowMethod
ddHodgeWorkflow(g, X, V;
    rdim=size(X,1)::Int, λ=0.1, ϵ=0.1,
    ssa=true, ssart=1, krytol=1e-32, useCUDA=false
) -> ddh::ddhResult

ddHodge standard workflow

Arguments

Required

  • g::SimpleGraph: undirected graph of N vertices.
  • X::AbstractMatrix: the matrices of data points.
  • V::AbstractMatrix: the matrices of velocites at the points.

These matrices X and V should be in the same size of M variables (rows) x N samples (columns) and the order of columns should correspond to the vertex order of g.

Optional

  • rdim::Int: the dimension of the tangent space. The default is M (no reduction).
  • λ::Float64=0.1: the regularization parameter for Hessian estimation.
  • ϵ::Float64=0.1: the regularization parameter for Jacobian estimation.
  • ssa::Bool=true: the flag of peforming subspace (axis) alignment.
  • ssart::Int=1: the reference node of axis alignment.
  • krytol::Float64=1e-32: the tolerance for Krylov solver.
  • useCUDA::Bool=false: the flag for GPU acceleration.

Output

The returned ddh::ddhResult consists of the following fields.

Vertex-level features

  • ddh.u: potential
  • ddh.div: divergence
  • ddh.rot: rotation (curl)
  • ddh.vgrass: Grassmann distance (averaged at vertices)
  • ddh.spans: local PCA tangent spaces
  • ddh.frames: axes-aligned local PCA spaces
  • ddh.planes: PC1-2 tangent planes
  • ddh.H: reconstructed Hessian matrices (in ddh.frames)
  • ddh.J: reconstructed Jacobian matrices (in ddh.frames)

Edge-level features

  • ddh.w: edge weight
  • ddh.egrass: Grassmann distance
  • ddh.rmaps: restriction maps

Operaters

  • ddh.d0: gradient operator
  • ddh.d0s: sheaf version of gradient operator

The order of these values is consistent with the vertex/edge order of given graph g.

source
ddHodge.graphgradMethod
graphgrad(g)

Construct grad operator from graph g

The direction is implicitly determined by node order, i.e., i -> j if i < j.

Example

julia> using ddHodge.Graphs

julia> g = complete_graph(3)
{3, 3} undirected simple Int64 graph

julia> collect(edges(g))
3-element Vector{Graphs.SimpleGraphs.SimpleEdge{Int64}}:
 Edge 1 => 2
 Edge 1 => 3
 Edge 2 => 3

julia> d0 = ddHodge.graphgrad(g)
3×3 SparseArrays.SparseMatrixCSC{Float64, Int64} with 6 stored entries:
 -1.0   1.0   ⋅ 
 -1.0    ⋅   1.0
   ⋅   -1.0  1.0
source
ddHodge.kNNGraphMethod
kNNGraph(X::Matrix, k::Int) -> (g::SimpleGraph, kdtree::KDTree)

Graph construction using k-NN

source
ddHodge.schurselectMethod
schurselect(F::Schur, idx::Vector{Int}) -> (values, Z::Matrix)

Returns the basis of invariant subspace Z; the columns are the selected Schur vectors specified by idx.

Warning

The invariance is guaranteed if the conjugate pair of the complex eigenvalues were selected.

Example

julia> using ddHodge.LinearAlgebra: schur

julia> A = [0 1 0; -1 0 0; 0 0 1]
3×3 Matrix{Int64}:
  0  1  0
 -1  0  0
  0  0  1

julia> F = schur(A);

julia> i, j, = sortperm(abs.(imag.(F.values)),rev=true);

julia> vals, Z = schurselect(F,[i,j]) # x-y plane is the invariant of A 
(values = ComplexF64[0.0 + 1.0im, 0.0 - 1.0im], Z = [1.0 0.0; 0.0 1.0; 0.0 0.0])
source