| Title: | Symmetric Matrix and Tensor Decomposition |
|---|---|
| Description: | Provides symmetric matrix and tensor operations and decomposition algorithms including symmetric NMF (symNMF), PageRank, Label Propagation, Higher-order Power Method, and TOPHITS. Designed to work with 'rTensor' objects. Methods are described in Kuang et al. (2012) <doi:10.1137/1.9781611972825.10>, Kolda and Mayo (2011) <doi:10.1137/100801482>, and Kolda and Bader (2006) <doi:10.1137/1.9781611972764.26>. |
| Authors: | Koki Tsuyuzaki [aut, cre] |
| Maintainer: | Koki Tsuyuzaki <[email protected]> |
| License: | MIT + file LICENSE |
| Version: | 0.1.0 |
| Built: | 2026-05-14 12:35:59 UTC |
| Source: | https://github.com/cran/symTensor |
Decomposes a symmetric tensor into a sum of rank-1 symmetric terms using the higher-order power method (S-HOPM) with deflation.
HigherOrderPower( x, R = 1L, max_iter = 1000L, tol = 1e-08, init = c("random", "svd") )HigherOrderPower( x, R = 1L, max_iter = 1000L, tol = 1e-08, init = c("random", "svd") )
x |
A symmetric rTensor Tensor object (all modes must be equal) |
R |
Number of rank-1 components to extract (default: 1) |
max_iter |
Maximum iterations per component (default: 1000) |
tol |
Convergence tolerance (default: 1e-8) |
init |
Initialization method: |
For a d-th order symmetric tensor ,
finds .
A list with components:
Numeric vector of length R with eigenvalues
Matrix (n x R) where each column is a unit eigenvector
Residual tensor after deflation
Logical vector of length R
Integer vector of iterations per component
Kolda, T. G., & Mayo, J. R. (2011). Shifted Power Method for Computing Tensor Eigenpairs. SIAM Journal on Matrix Analysis and Applications.
De Lathauwer, L., De Moor, B., & Vandewalle, J. (2000). On the Best Rank-1 and Rank-(R1, R2, ..., RN) Approximation of Higher-Order Tensors. SIAM Journal on Matrix Analysis and Applications.
library(rTensor) # Create a symmetric 3-mode tensor n <- 5 v1 <- rnorm(n); v1 <- v1 / sqrt(sum(v1^2)) # T = 3.0 * v1 o v1 o v1 arr <- 3.0 * outer(outer(v1, v1), v1) tnsr <- as.tensor(arr) result <- HigherOrderPower(tnsr, R = 1) result$lambdas # should be close to 3.0library(rTensor) # Create a symmetric 3-mode tensor n <- 5 v1 <- rnorm(n); v1 <- v1 / sqrt(sum(v1^2)) # T = 3.0 * v1 o v1 o v1 arr <- 3.0 * outer(outer(v1, v1), v1) tnsr <- as.tensor(arr) result <- HigherOrderPower(tnsr, R = 1) result$lambdas # should be close to 3.0
For a matrix, checks if A[i,j] == A[j,i] for all i,j. For a tensor, checks if the tensor is invariant under all permutations of its indices (i.e., fully symmetric).
isSymmetricTensor(x, tol = 1e-10)isSymmetricTensor(x, tol = 1e-10)
x |
A matrix or rTensor Tensor object |
tol |
Numeric tolerance for comparison (default: 1e-10) |
Logical indicating whether x is symmetric
# Matrix A <- matrix(c(1,2,2,3), 2, 2) isSymmetricTensor(A) # TRUE # 3-mode symmetric tensor library(rTensor) arr <- array(0, dim = c(3, 3, 3)) for (i in 1:3) for (j in 1:3) for (k in 1:3) { arr[i, j, k] <- i + j + k } isSymmetricTensor(as.tensor(arr)) # TRUE# Matrix A <- matrix(c(1,2,2,3), 2, 2) isSymmetricTensor(A) # TRUE # 3-mode symmetric tensor library(rTensor) arr <- array(0, dim = c(3, 3, 3)) for (i in 1:3) for (j in 1:3) for (k in 1:3) { arr[i, j, k] <- i + j + k } isSymmetricTensor(as.tensor(arr)) # TRUE
Semi-supervised label propagation on a graph. Uses Personalized PageRank internally to propagate labels from seed nodes to the rest of the graph.
LabelPropagation(A, seeds, damping = 0.85, max_iter = 1000L, tol = 1e-08)LabelPropagation(A, seeds, damping = 0.85, max_iter = 1000L, tol = 1e-08)
A |
Square non-negative adjacency/weight matrix (n x n) |
seeds |
Named integer vector or factor. Names correspond to node indices (1-based), values are the class labels. Unlabeled nodes will be classified. |
damping |
Damping factor passed to |
max_iter |
Maximum iterations for PageRank (default: 1000) |
tol |
Convergence tolerance for PageRank (default: 1e-8) |
A list with components:
Integer vector of predicted labels for all nodes
Matrix (n x K) of PageRank-based scores per class
Unique class labels
# 6-node graph with 2 communities A <- matrix(0, 6, 6) A[1,2] <- A[2,1] <- 1; A[1,3] <- A[3,1] <- 1; A[2,3] <- A[3,2] <- 1 A[4,5] <- A[5,4] <- 1; A[4,6] <- A[6,4] <- 1; A[5,6] <- A[6,5] <- 1 A[3,4] <- A[4,3] <- 0.5 # weak bridge seeds <- c("1" = 1, "6" = 2) # node 1 -> class 1, node 6 -> class 2 result <- LabelPropagation(A, seeds) result$labels# 6-node graph with 2 communities A <- matrix(0, 6, 6) A[1,2] <- A[2,1] <- 1; A[1,3] <- A[3,1] <- 1; A[2,3] <- A[3,2] <- 1 A[4,5] <- A[5,4] <- 1; A[4,6] <- A[6,4] <- 1; A[5,6] <- A[6,5] <- 1 A[3,4] <- A[4,3] <- 0.5 # weak bridge seeds <- c("1" = 1, "6" = 2) # node 1 -> class 1, node 6 -> class 2 result <- LabelPropagation(A, seeds) result$labels
Computes the PageRank vector for a given adjacency or transition matrix using power iteration.
PageRank( A, damping = 0.85, personalization = NULL, max_iter = 1000L, tol = 1e-08 )PageRank( A, damping = 0.85, personalization = NULL, max_iter = 1000L, tol = 1e-08 )
A |
Square non-negative matrix (adjacency or weight matrix). Columns are normalized internally to form a column-stochastic matrix. |
damping |
Damping factor (default: 0.85). Probability of following a link. |
personalization |
Optional personalization vector of length |
max_iter |
Maximum number of power iterations (default: 1000) |
tol |
Convergence tolerance on L1 norm change (default: 1e-8) |
A list with components:
PageRank vector (sums to 1)
Number of iterations until convergence
Logical, whether the algorithm converged
# Simple web graph A <- matrix(c(0,1,1,0, 1,0,0,1, 0,1,0,1, 1,0,1,0), 4, 4) pr <- PageRank(A) pr$vector# Simple web graph A <- matrix(c(0,1,1,0, 1,0,0,1, 0,1,0,1, 1,0,1,0), 4, 4) pr <- PageRank(A) pr$vector
For a matrix, computes (A + A^T) / 2. For a tensor, averages over all permutations of indices.
Symmetrize(x)Symmetrize(x)
x |
A matrix or rTensor Tensor object (must have equal modes) |
Symmetrized matrix or Tensor
A <- matrix(c(1, 2, 3, 4), 2, 2) Symmetrize(A) library(rTensor) arr <- array(rnorm(27), dim = c(3, 3, 3)) S <- Symmetrize(as.tensor(arr)) isSymmetricTensor(S) # TRUEA <- matrix(c(1, 2, 3, 4), 2, 2) Symmetrize(A) library(rTensor) arr <- array(rnorm(27), dim = c(3, 3, 3)) S <- Symmetrize(as.tensor(arr)) isSymmetricTensor(S) # TRUE
Decomposes a symmetric non-negative matrix
using multiplicative update (MU) rules with Beta-divergence.
symNMF( Omega, K, model = c("MSM", "MM"), Beta = 2, initM = NULL, initS = NULL, lambda = 0, num.iter = 100L, thr = 1e-10, verbose = FALSE )symNMF( Omega, K, model = c("MSM", "MM"), Beta = 2, initM = NULL, initS = NULL, lambda = 0, num.iter = 100L, thr = 1e-10, verbose = FALSE )
Omega |
Symmetric non-negative matrix (N x N) |
K |
Number of components (rank) |
model |
Character, either |
Beta |
Beta parameter for Beta-divergence (default: 2). 2 = Frobenius, 1 = KL, 0 = IS. |
initM |
Initial M matrix (N x K). If NULL, random non-negative initialization. |
initS |
Initial S matrix (K x K). If NULL, random non-negative initialization.
Ignored when |
lambda |
Regularization parameter for |
num.iter |
Maximum number of iterations (default: 100) |
thr |
Convergence threshold on relative change (default: 1e-10) |
verbose |
Logical, print iteration info (default: FALSE) |
Two models are supported:
MSM: where M is N x K and S is K x K
MM: where M is N x K (special case with S = I)
The loss function is the Beta-divergence :
: Frobenius (Euclidean) distance
: Kullback-Leibler divergence
: Itakura-Saito divergence
Other values of interpolate/extrapolate these cases
A list with components:
Factor matrix (N x K)
Coefficient matrix (K x K). Identity matrix when model = "MM".
Vector of reconstruction errors at each iteration
Vector of relative changes at each iteration
Logical
Number of iterations performed
Kuang, D., Park, H., & Ding, C. H. Q. (2012). Symmetric Nonnegative Matrix Factorization for Graph Clustering. SDM 2012.
Fevotte, C. & Idier, J. (2011). Algorithms for Nonnegative Matrix Factorization with the Beta-Divergence. Neural Computation, 23(9).
set.seed(42) M_true <- matrix(c(rep(c(1,0), 5), rep(c(0,1), 5)), 10, 2) Omega <- M_true %*% t(M_true) + matrix(runif(100, 0, 0.1), 10, 10) Omega <- (Omega + t(Omega)) / 2 # Frobenius (Beta=2, default) res_fro <- symNMF(Omega, K = 2, Beta = 2) # KL divergence (Beta=1) res_kl <- symNMF(Omega, K = 2, Beta = 1) # Itakura-Saito (Beta=0) res_is <- symNMF(Omega, K = 2, Beta = 0)set.seed(42) M_true <- matrix(c(rep(c(1,0), 5), rep(c(0,1), 5)), 10, 2) Omega <- M_true %*% t(M_true) + matrix(runif(100, 0, 0.1), 10, 10) Omega <- (Omega + t(Omega)) / 2 # Frobenius (Beta=2, default) res_fro <- symNMF(Omega, K = 2, Beta = 2) # KL divergence (Beta=1) res_kl <- symNMF(Omega, K = 2, Beta = 1) # Itakura-Saito (Beta=0) res_is <- symNMF(Omega, K = 2, Beta = 0)
Computes hub and authority scores for higher-order link analysis using Tucker decomposition of a tensor. Generalization of HITS (Hyperlink-Induced Topic Search) to tensors.
TOPHITS(x, R = 2L, max_iter = 100L, tol = 1e-08)TOPHITS(x, R = 2L, max_iter = 100L, tol = 1e-08)
x |
An rTensor Tensor object (2-mode or 3-mode) |
R |
Integer or integer vector specifying the number of components per mode. If scalar, same rank is used for all modes. |
max_iter |
Maximum iterations for Tucker decomposition (default: 100) |
tol |
Convergence tolerance (default: 1e-8) |
For a 3-mode tensor ,
TOPHITS computes Tucker decomposition
and derives hub/authority scores from the factor matrices.
For a symmetric tensor (all modes equal), a single set of scores is returned.
A list with components:
List of score matrices, one per mode (each n_i x R_i). For symmetric tensors, a single matrix.
Core tensor from Tucker decomposition
List of integer vectors giving node rankings per mode (ordered by dominant score)
Full Tucker decomposition result from rTensor
Kolda, T. G., & Bader, B. W. (2006). The TOPHITS Model for Higher-Order Web Link Analysis. Workshop on Link Analysis, Counterterrorism and Security.
library(rTensor) # Create a 3-mode tensor arr <- array(abs(rnorm(125)), dim = c(5, 5, 5)) tnsr <- as.tensor(arr) result <- TOPHITS(tnsr, R = 2) result$rankingslibrary(rTensor) # Create a 3-mode tensor arr <- array(abs(rnorm(125)), dim = c(5, 5, 5)) tnsr <- as.tensor(arr) result <- TOPHITS(tnsr, R = 2) result$rankings