Title: | Joint Sparse Optimization via Proximal Gradient Method for Cell Fate Conversion |
---|---|
Description: | Implementation of joint sparse optimization (JSparO) to infer the gene regulatory network for cell fate conversion. The proximal gradient method is implemented to solve different low-order regularization models for JSparO. |
Authors: | Xinlin Hu [aut, cre], Yaohua Hu [cph, aut] |
Maintainer: | Xinlin Hu <[email protected]> |
License: | GPL (>= 3) |
Version: | 1.5.0 |
Built: | 2024-12-01 07:59:43 UTC |
Source: | CRAN |
This is the main function of JSparO aimed to solve the low-order regularization models with norm.
demo_JSparO(A, B, X, s, p, q, maxIter = 200)
demo_JSparO(A, B, X, s, p, q, maxIter = 200)
A |
Gene expression data of transcriptome factors (i.e. feature matrix in machine learning). The dimension of A is m * n. |
B |
Gene expression data of target genes (i.e. observation matrix in machine learning). The dimension of B is m * t. |
X |
Gene expression data of Chromatin immunoprecipitation or other matrix (i.e. initial iterative point in machine learning). The dimension of X is n * t. |
s |
joint sparsity level |
p |
value for |
q |
value for |
maxIter |
maximum iteration |
The demo_JSparO function is used to solve joint sparse optimization problem via different algorithms.
Based on norm, functions with different p and q are implemented to solve the problem:
to obtain s-joint sparse solution.
The solution of proximal gradient method with regularizer.
Xinlin Hu [email protected]
Yaohua Hu [email protected]
m <- 256; n <- 1024; t <- 5; maxIter0 <- 50 A0 <- matrix(rnorm(m * n), nrow = m, ncol = n) B0 <- matrix(rnorm(m * t), nrow = m, ncol = t) X0 <- matrix(0, nrow = n, ncol = t) res_JSparO <- demo_JSparO(A0, B0, X0, s = 10, p = 2, q = 'half', maxIter = maxIter0)
m <- 256; n <- 1024; t <- 5; maxIter0 <- 50 A0 <- matrix(rnorm(m * n), nrow = m, ncol = n) B0 <- matrix(rnorm(m * t), nrow = m, ncol = t) X0 <- matrix(0, nrow = n, ncol = t) res_JSparO <- demo_JSparO(A0, B0, X0, s = 10, p = 2, q = 'half', maxIter = maxIter0)
normThe function aims to solve regularized least squares.
L1HalfThr(A, B, X, s, maxIter = 200)
L1HalfThr(A, B, X, s, maxIter = 200)
A |
Gene expression data of transcriptome factors (i.e. feature matrix in machine learning). The dimension of A is m * n. |
B |
Gene expression data of target genes (i.e. observation matrix in machine learning). The dimension of B is m * t. |
X |
Gene expression data of Chromatin immunoprecipitation or other matrix (i.e. initial iterative point in machine learning). The dimension of X is n * t. |
s |
joint sparsity level |
maxIter |
maximum iteration |
The L1HalfThr function aims to solve the problem:
to obtain s-joint sparse solution.
The solution of proximal gradient method with regularizer.
Xinlin Hu [email protected]
Yaohua Hu [email protected]
m <- 256; n <- 1024; t <- 5; maxIter0 <- 50 A0 <- matrix(rnorm(m * n), nrow = m, ncol = n) B0 <- matrix(rnorm(m * t), nrow = m, ncol = t) X0 <- matrix(0, nrow = n, ncol = t) NoA <- norm(A0, '2'); A0 <- A0/NoA; B0 <- B0/NoA res_L1half <- L1HalfThr(A0, B0, X0, s = 10, maxIter = maxIter0)
m <- 256; n <- 1024; t <- 5; maxIter0 <- 50 A0 <- matrix(rnorm(m * n), nrow = m, ncol = n) B0 <- matrix(rnorm(m * t), nrow = m, ncol = t) X0 <- matrix(0, nrow = n, ncol = t) NoA <- norm(A0, '2'); A0 <- A0/NoA; B0 <- B0/NoA res_L1half <- L1HalfThr(A0, B0, X0, s = 10, maxIter = maxIter0)
normThe function aims to solve regularized least squares.
L1HardThr(A, B, X, s, maxIter = 200)
L1HardThr(A, B, X, s, maxIter = 200)
A |
Gene expression data of transcriptome factors (i.e. feature matrix in machine learning). The dimension of A is m * n. |
B |
Gene expression data of target genes (i.e. observation matrix in machine learning). The dimension of B is m * t. |
X |
Gene expression data of Chromatin immunoprecipitation or other matrix (i.e. initial iterative point in machine learning). The dimension of X is n * t. |
s |
joint sparsity level |
maxIter |
maximum iteration |
The L1HardThr function aims to solve the problem:
to obtain s-joint sparse solution.
The solution of proximal gradient method with regularizer.
Xinlin Hu [email protected]
Yaohua Hu [email protected]
m <- 256; n <- 1024; t <- 5; maxIter0 <- 50 A0 <- matrix(rnorm(m * n), nrow = m, ncol = n) B0 <- matrix(rnorm(m * t), nrow = m, ncol = t) X0 <- matrix(0, nrow = n, ncol = t) NoA <- norm(A0, '2'); A0 <- A0/NoA; B0 <- B0/NoA res_L10 <- L1HardThr(A0, B0, X0, s = 10, maxIter = maxIter0)
m <- 256; n <- 1024; t <- 5; maxIter0 <- 50 A0 <- matrix(rnorm(m * n), nrow = m, ncol = n) B0 <- matrix(rnorm(m * t), nrow = m, ncol = t) X0 <- matrix(0, nrow = n, ncol = t) NoA <- norm(A0, '2'); A0 <- A0/NoA; B0 <- B0/NoA res_L10 <- L1HardThr(A0, B0, X0, s = 10, maxIter = maxIter0)
The function aims to compute the norm.
L1normFun(x)
L1normFun(x)
x |
vector |
The L1normFun aims to compute the norm:
The norm of vector x
Xinlin Hu [email protected]
Yaohua Hu [email protected]
normThe function aims to solve regularized least squares.
L1SoftThr(A, B, X, s, maxIter = 200)
L1SoftThr(A, B, X, s, maxIter = 200)
A |
Gene expression data of transcriptome factors (i.e. feature matrix in machine learning). The dimension of A is m * n. |
B |
Gene expression data of target genes (i.e. observation matrix in machine learning). The dimension of B is m * t. |
X |
Gene expression data of Chromatin immunoprecipitation or other matrix (i.e. initial iterative point in machine learning). The dimension of X is n * t. |
s |
joint sparsity level |
maxIter |
maximum iteration |
The L1SoftThr function aims to solve the problem:
to obtain s-joint sparse solution.
The solution of proximal gradient method with regularizer.
Xinlin Hu [email protected]
Yaohua Hu [email protected]
m <- 256; n <- 1024; t <- 5; maxIter0 <- 50 A0 <- matrix(rnorm(m * n), nrow = m, ncol = n) B0 <- matrix(rnorm(m * t), nrow = m, ncol = t) X0 <- matrix(0, nrow = n, ncol = t) NoA <- norm(A0, '2'); A0 <- A0/NoA; B0 <- B0/NoA res_L11 <- L1SoftThr(A0, B0, X0, s = 10, maxIter = maxIter0)
m <- 256; n <- 1024; t <- 5; maxIter0 <- 50 A0 <- matrix(rnorm(m * n), nrow = m, ncol = n) B0 <- matrix(rnorm(m * t), nrow = m, ncol = t) X0 <- matrix(0, nrow = n, ncol = t) NoA <- norm(A0, '2'); A0 <- A0/NoA; B0 <- B0/NoA res_L11 <- L1SoftThr(A0, B0, X0, s = 10, maxIter = maxIter0)
normThe function aims to solve regularized least squares.
L1twothirdsThr(A, B, X, s, maxIter = 200)
L1twothirdsThr(A, B, X, s, maxIter = 200)
A |
Gene expression data of transcriptome factors (i.e. feature matrix in machine learning). The dimension of A is m * n. |
B |
Gene expression data of target genes (i.e. observation matrix in machine learning). The dimension of B is m * t. |
X |
Gene expression data of Chromatin immunoprecipitation or other matrix (i.e. initial iterative point in machine learning). The dimension of X is n * t. |
s |
joint sparsity level |
maxIter |
maximum iteration |
The L1twothirdsThr function aims to solve the problem:
to obtain s-joint sparse solution.
The solution of proximal gradient method with regularizer.
Xinlin Hu [email protected]
Yaohua Hu [email protected]
m <- 256; n <- 1024; t <- 5; maxIter0 <- 50 A0 <- matrix(rnorm(m * n), nrow = m, ncol = n) B0 <- matrix(rnorm(m * t), nrow = m, ncol = t) X0 <- matrix(0, nrow = n, ncol = t) NoA <- norm(A0, '2'); A0 <- A0/NoA; B0 <- B0/NoA res_L1twothirds <- L1twothirdsThr(A0, B0, X0, s = 10, maxIter = maxIter0)
m <- 256; n <- 1024; t <- 5; maxIter0 <- 50 A0 <- matrix(rnorm(m * n), nrow = m, ncol = n) B0 <- matrix(rnorm(m * t), nrow = m, ncol = t) X0 <- matrix(0, nrow = n, ncol = t) NoA <- norm(A0, '2'); A0 <- A0/NoA; B0 <- B0/NoA res_L1twothirds <- L1twothirdsThr(A0, B0, X0, s = 10, maxIter = maxIter0)
normThe function aims to solve regularized least squares.
L2HalfThr(A, B, X, s, maxIter = 200)
L2HalfThr(A, B, X, s, maxIter = 200)
A |
Gene expression data of transcriptome factors (i.e. feature matrix in machine learning). The dimension of A is m * n. |
B |
Gene expression data of target genes (i.e. observation matrix in machine learning). The dimension of B is m * t. |
X |
Gene expression data of Chromatin immunoprecipitation or other matrix (i.e. initial iterative point in machine learning). The dimension of X is n * t. |
s |
joint sparsity level |
maxIter |
maximum iteration |
The L2HalfThr function aims to solve the problem:
to obtain s-joint sparse solution.
The solution of proximal gradient method with regularizer.
Xinlin Hu [email protected]
Yaohua Hu [email protected]
m <- 256; n <- 1024; t <- 5; maxIter0 <- 50 A0 <- matrix(rnorm(m * n), nrow = m, ncol = n) B0 <- matrix(rnorm(m * t), nrow = m, ncol = t) X0 <- matrix(0, nrow = n, ncol = t) NoA <- norm(A0, '2'); A0 <- A0/NoA; B0 <- B0/NoA res_L2half <- L2HalfThr(A0, B0, X0, s = 10, maxIter = maxIter0)
m <- 256; n <- 1024; t <- 5; maxIter0 <- 50 A0 <- matrix(rnorm(m * n), nrow = m, ncol = n) B0 <- matrix(rnorm(m * t), nrow = m, ncol = t) X0 <- matrix(0, nrow = n, ncol = t) NoA <- norm(A0, '2'); A0 <- A0/NoA; B0 <- B0/NoA res_L2half <- L2HalfThr(A0, B0, X0, s = 10, maxIter = maxIter0)
normThe function aims to solve regularized least squares.
L2HardThr(A, B, X, s, maxIter = 200)
L2HardThr(A, B, X, s, maxIter = 200)
A |
Gene expression data of transcriptome factors (i.e. feature matrix in machine learning). The dimension of A is m * n. |
B |
Gene expression data of target genes (i.e. observation matrix in machine learning). The dimension of B is m * t. |
X |
Gene expression data of Chromatin immunoprecipitation or other matrix (i.e. initial iterative point in machine learning). The dimension of X is n * t. |
s |
joint sparsity level |
maxIter |
maximum iteration |
The L2HardThr function aims to solve the problem:
to obtain s-joint sparse solution.
The solution of proximal gradient method with regularizer.
Xinlin Hu [email protected]
Yaohua Hu [email protected]
m <- 256; n <- 1024; t <- 5; maxIter0 <- 50 A0 <- matrix(rnorm(m * n), nrow = m, ncol = n) B0 <- matrix(rnorm(m * t), nrow = m, ncol = t) X0 <- matrix(0, nrow = n, ncol = t) NoA <- norm(A0, '2'); A0 <- A0/NoA; B0 <- B0/NoA res_L20 <- L2HardThr(A0, B0, X0, s = 10, maxIter = maxIter0)
m <- 256; n <- 1024; t <- 5; maxIter0 <- 50 A0 <- matrix(rnorm(m * n), nrow = m, ncol = n) B0 <- matrix(rnorm(m * t), nrow = m, ncol = t) X0 <- matrix(0, nrow = n, ncol = t) NoA <- norm(A0, '2'); A0 <- A0/NoA; B0 <- B0/NoA res_L20 <- L2HardThr(A0, B0, X0, s = 10, maxIter = maxIter0)
norm with Newton methodThe function aims to solve regularized least squares, where the proximal optimization subproblems will be solved by Newton method.
L2NewtonThr(A, B, X, s, q, maxIter = 200, innMaxIter = 30, innEps = 1e-06)
L2NewtonThr(A, B, X, s, q, maxIter = 200, innMaxIter = 30, innEps = 1e-06)
A |
Gene expression data of transcriptome factors (i.e. feature matrix in machine learning). The dimension of A is m * n. |
B |
Gene expression data of target genes (i.e. observation matrix in machine learning). The dimension of B is m * t. |
X |
Gene expression data of Chromatin immunoprecipitation or other matrix (i.e. initial iterative point in machine learning). The dimension of X is n * t. |
s |
joint sparsity level |
q |
value for |
maxIter |
maximum iteration |
innMaxIter |
maximum iteration in Newton step |
innEps |
criterion to stop inner iteration |
The L2NewtonThr function aims to solve the problem:
to obtain s-joint sparse solution.
The solution of proximal gradient method with regularizer.
Xinlin Hu [email protected]
Yaohua Hu [email protected]
m <- 256; n <- 1024; t <- 5; maxIter0 <- 50 A0 <- matrix(rnorm(m * n), nrow = m, ncol = n) B0 <- matrix(rnorm(m * t), nrow = m, ncol = t) X0 <- matrix(0, nrow = n, ncol = t) NoA <- norm(A0, '2'); A0 <- A0/NoA; B0 <- B0/NoA res_L2q <- L2NewtonThr(A0, B0, X0, s = 10, q = 0.2, maxIter = maxIter0)
m <- 256; n <- 1024; t <- 5; maxIter0 <- 50 A0 <- matrix(rnorm(m * n), nrow = m, ncol = n) B0 <- matrix(rnorm(m * t), nrow = m, ncol = t) X0 <- matrix(0, nrow = n, ncol = t) NoA <- norm(A0, '2'); A0 <- A0/NoA; B0 <- B0/NoA res_L2q <- L2NewtonThr(A0, B0, X0, s = 10, q = 0.2, maxIter = maxIter0)
normThe function aims to solve regularized least squares.
L2SoftThr(A, B, X, s, maxIter = 200)
L2SoftThr(A, B, X, s, maxIter = 200)
A |
Gene expression data of transcriptome factors (i.e. feature matrix in machine learning). The dimension of A is m * n. |
B |
Gene expression data of target genes (i.e. observation matrix in machine learning). The dimension of B is m * t. |
X |
Gene expression data of Chromatin immunoprecipitation or other matrix (i.e. initial iterative point in machine learning). The dimension of X is n * t. |
s |
joint sparsity level |
maxIter |
maximum iteration |
The L2SoftThr function aims to solve the problem:
to obtain s-joint sparse solution.
The solution of proximal gradient method with regularizer.
Xinlin Hu [email protected]
Yaohua Hu [email protected]
m <- 256; n <- 1024; t <- 5; maxIter0 <- 50 A0 <- matrix(rnorm(m * n), nrow = m, ncol = n) B0 <- matrix(rnorm(m * t), nrow = m, ncol = t) X0 <- matrix(0, nrow = n, ncol = t) NoA <- norm(A0, '2'); A0 <- A0/NoA; B0 <- B0/NoA res_L21 <- L2SoftThr(A0, B0, X0, s = 10, maxIter = maxIter0)
m <- 256; n <- 1024; t <- 5; maxIter0 <- 50 A0 <- matrix(rnorm(m * n), nrow = m, ncol = n) B0 <- matrix(rnorm(m * t), nrow = m, ncol = t) X0 <- matrix(0, nrow = n, ncol = t) NoA <- norm(A0, '2'); A0 <- A0/NoA; B0 <- B0/NoA res_L21 <- L2SoftThr(A0, B0, X0, s = 10, maxIter = maxIter0)
normThe function aims to solve regularized least squares.
L2twothirdsThr(A, B, X, s, maxIter = 200)
L2twothirdsThr(A, B, X, s, maxIter = 200)
A |
Gene expression data of transcriptome factors (i.e. feature matrix in machine learning). The dimension of A is m * n. |
B |
Gene expression data of target genes (i.e. observation matrix in machine learning). The dimension of B is m * t. |
X |
Gene expression data of Chromatin immunoprecipitation or other matrix (i.e. initial iterative point in machine learning). The dimension of X is n * t. |
s |
joint sparsity level |
maxIter |
maximum iteration |
The L2twothirdsThr function aims to solve the problem:
to obtain s-joint sparse solution.
The solution of proximal gradient method with regularizer.
Xinlin Hu [email protected]
Yaohua Hu [email protected]
m <- 256; n <- 1024; t <- 5; maxIter0 <- 50 A0 <- matrix(rnorm(m * n), nrow = m, ncol = n) B0 <- matrix(rnorm(m * t), nrow = m, ncol = t) X0 <- matrix(0, nrow = n, ncol = t) NoA <- norm(A0, '2'); A0 <- A0/NoA; B0 <- B0/NoA res_L2twothirds <- L2twothirdsThr(A0, B0, X0, s = 10, maxIter = maxIter0)
m <- 256; n <- 1024; t <- 5; maxIter0 <- 50 A0 <- matrix(rnorm(m * n), nrow = m, ncol = n) B0 <- matrix(rnorm(m * t), nrow = m, ncol = t) X0 <- matrix(0, nrow = n, ncol = t) NoA <- norm(A0, '2'); A0 <- A0/NoA; B0 <- B0/NoA res_L2twothirds <- L2twothirdsThr(A0, B0, X0, s = 10, maxIter = maxIter0)