| Title: | Inferring Functional Gene Co-Expression Networks from Single Cell Data |
|---|---|
| Description: | Uses statistical network modeling to understand the co-expression relationships among genes and to construct sparse gene co-expression networks from single-cell gene expression data. |
| Authors: | Wei Vivian Li [aut, cre]
|
| Maintainer: | Wei Vivian Li <[email protected]> |
| License: | GPL-3 |
| Version: | 1.0.1 |
| Built: | 2026-05-10 05:37:52 UTC |
| Source: | https://github.com/cran/scLink |
Calculate scLink's correlation matrix
sclink_cor(expr, ncores, nthre = 20, dthre = 0.9)sclink_cor(expr, ncores, nthre = 20, dthre = 0.9)
expr |
A gene expression matrix with rows representing cells and columns representing genes.
Gene names are given as column names. Can be the output of |
ncores |
Number of cores if using parallel computation. |
nthre |
An integer specifying a threshold on the number of complete observations. Defaults to 20. |
dthre |
A number specifying the threshold on dropout probabilities. Defaults to 0.9. |
A correlation matrix for gene co-expression relationships.
Wei Vivian Li, [email protected]
count = readRDS(system.file("extdata", "example.rds", package = "scLink")) count.norm = sclink_norm(count, scale.factor = 1e6, filter.genes = TRUE, n = 500) corr = sclink_cor(expr = count.norm, ncores = 1)count = readRDS(system.file("extdata", "example.rds", package = "scLink")) count.norm = sclink_norm(count, scale.factor = 1e6, filter.genes = TRUE, n = 500) corr = sclink_cor(expr = count.norm, ncores = 1)
Infer gene co-expression networks
sclink_net(expr, ncores, lda = seq(1, 0.1, -0.05), nthre = 20, dthre = 0.9)sclink_net(expr, ncores, lda = seq(1, 0.1, -0.05), nthre = 20, dthre = 0.9)
expr |
A gene expression matrix with rows representing cells and columns representing genes.
Gene names are given as column names. Can be the output of |
ncores |
Number of cores if using parallel computation. |
lda |
A vector specifying a sequence of lambda values to be used in the penalized likelihood. |
nthre |
An integer specifying a threshold on the number of complete observations. Defaults to 20. |
dthre |
A number specifying the threshold on dropout probabilities. Defaults to 0.9. |
A list for gene co-expression relationships. The list contains a cor element for
scLink's correlation matrix and a summary element for the gene networks. summary is a list
with each element corresponding to the result of one lambda value. Each element of summary
contains the following information:
the adjacency matrix specifying the gene-gene edges;
the estimated concentration matrix;
number of edges in the gene network;
BIC score;
value of lambda in the penalty.
Wei Vivian Li, [email protected]
count = readRDS(system.file("extdata", "example.rds", package = "scLink")) count.norm = sclink_norm(count, scale.factor = 1e6, filter.genes = TRUE, n = 500) networks = sclink_net(expr = count.norm, ncores = 1, lda = seq(0.5, 0.1, -0.05))count = readRDS(system.file("extdata", "example.rds", package = "scLink")) count.norm = sclink_norm(count, scale.factor = 1e6, filter.genes = TRUE, n = 500) networks = sclink_net(expr = count.norm, ncores = 1, lda = seq(0.5, 0.1, -0.05))
Pre-process data for scLink
sclink_norm( count, scale.factor = 1e+06, filter.genes = FALSE, gene.names = NULL, n = 500 )sclink_norm( count, scale.factor = 1e+06, filter.genes = FALSE, gene.names = NULL, n = 500 )
count |
A full gene count matrix with rows representing cells and columns representing genes. Gene names are given as column names. |
scale.factor |
A number specifying the sclae factor used for library size normalization. Defaults to 1e6. |
filter.genes |
A Boolean specifying whether scLink should select genes based on mean expression.
When set to |
gene.names |
A character vector specifying the genes used for network construction.
Only needed when |
n |
An integer specifying the number of genes to be selected by scLink (defaults to 500).
Only needed when |
A transformed and normalized gene expression matrix that can be used for correlation calculation and network construction.
Wei Vivian Li, [email protected]
count = readRDS(system.file("extdata", "example.rds", package = "scLink")) count.norm = sclink_norm(count, scale.factor = 1e6, filter.genes = TRUE, n = 500)count = readRDS(system.file("extdata", "example.rds", package = "scLink")) count.norm = sclink_norm(count, scale.factor = 1e6, filter.genes = TRUE, n = 500)