Index

TemporalNetworks.BlockGraphType
BlockGraph(N :: Int64, T :: Int64, list :: Vector, η :: Float64, clusters :: Union{Nothing, Vector}, degrees :: Union{Nothing, Vector})

Multiplex type block graph instance. Constructs a temporal network of N vertices per time slice with T time steps.

The transition types are encoded in list which is a Vector{Int64}. Numbers indicate number of disjoint clusters to be constructed at certain time points. clusters and degrees hold information about which vertices belong to which clusters and the cluster degrees respectively (the clusters are regular subgraphs).

The instance used as a function and passed with no arguments returns a sequence of adjacency matrices of type Vector{Float64}.

source
TemporalNetworks.BlockGraphNonMultiplexType
BlockGraphNonMultiplex(N :: Int64, T :: Int64, list :: Vector, η :: Float64, clusters :: Union{Nothing, Vector}, degrees :: Union{Nothing, Vector}, interaction :: IntraBlockConnectivity, evolve :: Int64)

Non-multiplex type staircase graph instance. Constructs a temporal network of Nspatial vertices $\bigcup_{t=1}^T\{x: (t,x)\in\mathcal{V}\}$ with T time steps.

The transition types are encoded in list which is a Vector{Int64}. Numbers indicate number of disjoint clusters to be constructed at certain time points. clusters and degrees hold information about which vertices belong to which clusters and the cluster degrees respectively (the clusters are regular subgraphs).

interaction is a IntraBlockConnectivity instance that acts as a function, and is applied to all intercluster edge weights. The default is ScalingConnectivity(T) which is the function $x \mapsto x/T$.

evolve indicates how many vertices per block are to be considered absent in the next layer, and switched with equally many new vertices from the absent ones. Defaults to 1.

BlockGraphNonMultiplex(N, T, list, η, clusters, degrees)

BlockGraphNonMultiplex instance with interactivity = ScalingConnectivity(2) and evolve=1.

BlockGraphNonMultiplex(N, T, list, η, clusters, degrees)

BlockGraphNonMultiplex instance with interactivity = ScalingConnectivity(2)

The instance used as a function and passed with no arguments returns a sequence of adjacency matrices of type Vector{Float64}.

source
TemporalNetworks.DegreeNormalizationType
(:: DegreeNormalization)(L :: Matrix)

Normalizes Laplacian with resepct to the standard normalization $D^{-1/2}LD^{-1/2}$ where D is diagonal with i-th entry equal to deg(i).

source
TemporalNetworks.MultilayerGraphType
MultilayerGraph{M <: TemporalConnectivity}(W :: Vector{Matrix{Float64}}; connect = Multiplex())

The temporal network instance. Stores sequence of adjacencies W, nodes per layer N, number of layers T and the connectivity of abstract type TemporalConnectivity. Default is Multiplex(), other options include NonMultiplexCompressed().

source
TemporalNetworks.MultiplexType
(:: Multiplex)(mlgraph :: MultilayerGraph)

Connection type :: TemporalConnectivity for MultilayerGraph. Returns W_spat, W_temp such that $Wˢᵖᵃᵗ = ⨁ᵢWⁱ$ and $Wᵗᵉᵐᵖ = W'⊗Iₙ$.

Equivalently,

  1. W_spat = directsum(mlgraph.W) where W :: Vector{Matrix}
  2. W_temp = kron(Wt, Matrix{Float64}(I, mlgraph.N, mlgraph.N)) where Wt is a T x T matrix describing how the layers are connected (symmetric supradiagonal).
source
TemporalNetworks.NonMultiplexCompressedType
(:: NonMultiplexCompressed)(mlgraph :: MultilayerGraph)

Nonmultiplex connection type :: TemporalConnectivity for MultilayerGraph. Returns W_spat, W_temp such that:

  1. W_spat = directsum(mlgraph.W) where W :: Vector{Matrix}
  2. W_temp[i,j] ≂̸ 0 iff
    • mlgraph.W[k][q,:] ≂̸ 0
    • mlgraph.W[l][r,:] ≂̸ 0
    where $(k,q) ∈ Rᴺ × Rᵀ$ and $(l,r) ∈ Rᴺ × Rᵀ$ are the separated spacetime indices corresponding to $i ∈ Rᴺᵀ$ and $j ∈ Rᴺᵀ$ respectively.
source
TemporalNetworks.RayleighBalancingType
RayleighBalancing(ind :: Int64)

Creates DiffusionEstimator instance to compute the diffusion constant $a$ for the supra-Laplacian arising from a non-multiplex network by matching the spatial and temproral Rayleigh quotient contributions with respect to the eigenvalue indexed by ind.

Recommendeded for nonmultiplex networks.

(rb:: RayleighBalancing)(mlgraph :: MultilayerGraph ,norm :: Normalization, L_spat :: Matrix{Float64}, L_temp :: Matrix{Float64})

Computes the diffusion constant $a$ using the Rayleigh balancing criterion. See [FroylandKaliaKoltai2024] for details on this heuristic.

source
TemporalNetworks.SEBAPartitionType
SEBAPartition{MG <: MultilayerGraph, N <: Normalization}(partition :: SpectralPartition, )

SEBA partitiom instance. Stores the SpectralPartition instance in partition, supra-Laplacian eigenvector indices used in SEBA in inds, SEBA vectors in vecs and corresponding Cheeger ratios in cuts (depending on type of Normalization).

SEBAPartition applies the SEBA algorithm to chosen/identified eigenvectors from partition.evecs and computes Cheeger rations.

SEBAPartition(partition :: SpectralPartition, inds :: Vector{Int64})

Computes SEBAPartition instance using partition.evecs[:, inds].

SEBAPartition(partition :: SpectralPartition, max_ind :: Int64)

Computes SEBAPartition instance by finding max_ind leading non-trivial spatial eigenvectors within partition.evecs. Does not work with NonMultiplex type graphs.

source
TemporalNetworks.SpatTempMatchingType
(:: SpatTempMatching)(::SpatTempMatching)(mlgraph :: MultilayerGraph{T}, norm :: Normalization, L_spat :: Matrix{Float64}, L_temp :: Matrix{Float64}) where T <: Multiplex

Computes the diffusion constant $a$ by matching the second spatial eigenvalue (first non-trivial spatial eigenvalue) to the first temporal eigenvalue for the supra-Laplacian corresponding to a multiplex temporal network.

The multiplex structure guarantees a unique smallest intersection point. See [FroylandKoltai2023] and [AtnipFroylandKoltai2024] for details on this heuristic.

source
TemporalNetworks.SpatTempMatchingNonMultiplexType
SpatTempMatchingNonMultiplex(evec_ind :: Int64, temp_ind :: Union{Int64, Nothng}, thresh :: Float64)

Creates DiffusionEstimator instance to compute the diffusion constant $a$ for the supra-Laplacian arising from a non-multiplex network by matching the eigenvalue corresponding to evec_ind to the corresponding multiplex temporal eigenvalue with index temp_ind.

The parameter thresh is used by isspat(::MultilayerGraph{NonMultiplex}, ...) to check for spatial-like behaviour.

(:: SpatTempMatchingNonMultiplex)(::SpatTempMatching)(mlgraph :: MultilayerGraph{T}, norm :: Normalization, L_spat :: Matrix{Float64}, L_temp :: Matrix{Float64}) where T <: NonMultiplex

Computes the diffusion constant $a$ using the matching heuristic for non-multiplex networks. See [FroylandKaliaKoltai2024] for details on this heuristic.

source
TemporalNetworks.SpectralPartitionType
SpectralPartition{MG <: MultilayerGraph, N <: Normalization}(mlgraph :: MultilayerGraph; compute_a = SpatTempMatching(), norm = IdentityNormalization())

Spectral partition instance for the temporal network defined by mlgraph :: MultilayerGraph. mlgraph is stored in graph, norm :: Normalization stores the Laplacian normalization function (default is IdentityNormalization()), a contains the diffusion constant, L_temp and L_spat are the temporal and spatial supra-Laplacians respectively and evecs and evals store the eigendata of the supra-Laplacian.

SpectralPartition computes the supra-Laplacian $\mathbf L^{(a)} = \mathbf L^{\rm spat} + a^2 \mathbf L^{\rm temp}$, computes the appropriate value of the diffusion constant $a$ and computes the spectrum for both multiplex and nonmultiplex type networks.

Note that for nonmultiplex networks, evecs are lifted to $\mathbb R^{TN}$.

source
TemporalNetworks.bipartite_minimum_weight_edge_coverMethod

Computes the minimum-weight edge cover of a bipartite graph G. We assume there are m nodes on the left and n nodes on the right, where m≥n. The (i,j)th entry of the m×n matrix G contains the weight of the edge joining node i and node j. A binary m×n array is output, encoding the cover.

source
TemporalNetworks.bisectionMethod
bisection(f, a, b; fa = f(a), fb = f(b), ftol, wtol)

Bisection algorithm for finding the root $f(x) ≈ 0$ within the initial bracket [a,b].

Returns a named tuple (x = x, fx = f(x), isroot = ::Bool, iter = ::Int, ismaxiter = ::Bool).

Terminates when either

  1. abs(f(x)) < ftol (isroot = true),
  2. the width of the bracket is ≤wtol (isroot = true, to account for discontinuities),
  3. maxiter number of iterations is reached. (isroot = false, maxiter = true).

which are tested for in the above order. Therefore, care should be taken not to make wtol too large.

source