Title: | scRecover for imputation of single-cell RNA-seq data |
---|---|
Description: | scRecover is an R package for imputation of single-cell RNA-seq (scRNA-seq) data. It will detect and impute dropout values in a scRNA-seq raw read counts matrix while keeping the real zeros unchanged, since there are both dropout zeros and real zeros in scRNA-seq data. By combination with scImpute, SAVER and MAGIC, scRecover not only detects dropout and real zeros at higher accuracy, but also improve the downstream clustering and visualization results. |
Authors: | Zhun Miao, Xuegong Zhang <[email protected]> |
Maintainer: | Zhun Miao <[email protected]> |
License: | GPL |
Version: | 1.23.0 |
Built: | 2024-11-07 06:15:09 UTC |
Source: | https://github.com/bioc/scRecover |
This function is used to downsample the read counts in a cell for single-cell RNA-seq (scRNA-seq) data. It takes a non-negative vector of scRNA-seq raw read counts of a cell as input.
countsSampling(counts, fraction = 0.1)
countsSampling(counts, fraction = 0.1)
counts |
A cell's raw read counts for each gene, could be a vector or a SingleCellExperiment object. |
fraction |
Fraction of reads to be downsampled, should be between 0-1, default is 0.1. |
A vector of the downsampled read counts of each gene in the cell.
Zhun Miao.
scRecover
, for imputation of single-cell RNA-seq data.
estDropoutNum
, for estimating dropout gene number in a cell.
normalization
, for normalization of single-cell RNA-seq data.
scRecoverTest
, a test dataset for scRecover.
# Load test data data(scRecoverTest) # Downsample the read counts in oneCell oneCell.down <- countsSampling(counts = oneCell, fraction = 0.1)
# Load test data data(scRecoverTest) # Downsample the read counts in oneCell oneCell.down <- countsSampling(counts = oneCell, fraction = 0.1)
This function is used to estimate dropout gene number in a cell for single-cell RNA-seq (scRNA-seq) data. It takes a non-negative vector of scRNA-seq raw read counts of a cell as input.
estDropoutNum(sample = NULL, depth = 20, histCounts = NULL, return = "dropoutNum")
estDropoutNum(sample = NULL, depth = 20, histCounts = NULL, return = "dropoutNum")
sample |
A cell's raw read counts for each gene, could be a vector or a SingleCellExperiment object. |
depth |
Relative sequencing depth to be predicted compared with initial sample depth, should between 0-100, default is 20. |
histCounts |
Optional. Only needed when |
return |
A character for choosing the return value type of the function. "dropoutNum" (default) for dropout gene number, "geneNumPredict" for all expressed gene number predicted, "transcriptNum" for all transcript number predicted. |
The dropout gene number (or all expressed gene number) predicted in a cell.
Zhun Miao.
scRecover
, for imputation of single-cell RNA-seq data.
countsSampling
, for downsampling the read counts in a cell.
normalization
, for normalization of single-cell RNA-seq data.
scRecoverTest
, a test dataset for scRecover.
# Load test data data(scRecoverTest) # Estimate dropout gene number in a cell estDropoutNum(sample = counts[,1], return = "dropoutNum") # Estimate all expressed gene number in a cell estDropoutNum(sample = counts[,1], return = "geneNumPredict")
# Load test data data(scRecoverTest) # Estimate dropout gene number in a cell estDropoutNum(sample = counts[,1], return = "dropoutNum") # Estimate all expressed gene number in a cell estDropoutNum(sample = counts[,1], return = "geneNumPredict")
This function is used to normalize single-cell RNA-seq (scRNA-seq) data. It takes a non-negative matrix of scRNA-seq raw read counts or a SingleCellExperiment
object as input.
normalization(counts)
normalization(counts)
counts |
A non-negative integer matrix of scRNA-seq raw read counts or a |
A normalized scRNA-seq read counts matrix.
Zhun Miao.
scRecover
, for imputation of single-cell RNA-seq data.
estDropoutNum
, for estimating dropout gene number in a cell.
countsSampling
, for downsampling the read counts in a cell.
scRecoverTest
, a test dataset for scRecover.
# Load test data data(scRecoverTest) # Normalization of counts counts.norm <- normalization(counts = counts)
# Load test data data(scRecoverTest) # Normalization of counts counts.norm <- normalization(counts = counts)
This function is used to impute missing values in single-cell RNA-seq (scRNA-seq) data. It takes a non-negative matrix of scRNA-seq raw read counts or a SingleCellExperiment
object as input. So users should map the reads (obtained from sequencing libraries of the samples) to the corresponding genome and count the reads mapped to each gene according to the gene annotation to get the raw read counts matrix in advance.
scRecover(counts, Kcluster = NULL, labels = NULL, outputDir = NULL, depth = 20, SAVER = FALSE, MAGIC = FALSE, UMI = FALSE, hist_raw_counts = NULL, hist_RUG_counts = NULL, parallel = FALSE, BPPARAM = bpparam(), verbose = TRUE)
scRecover(counts, Kcluster = NULL, labels = NULL, outputDir = NULL, depth = 20, SAVER = FALSE, MAGIC = FALSE, UMI = FALSE, hist_raw_counts = NULL, hist_RUG_counts = NULL, parallel = FALSE, BPPARAM = bpparam(), verbose = TRUE)
counts |
A non-negative integer matrix of scRNA-seq raw read counts or a |
Kcluster |
An integer specifying the number of cell subpopulations. This parameter can be determined based on prior knowledge or clustering of raw data. |
labels |
Optional. Only needed when |
outputDir |
The path of the output directory. If not specified, a folder named with prefix 'outDir_scRecover_' under the temporary directory will be used. |
depth |
Relative sequencing depth to be predicted compared with initial sample depth, should between 2-100, default is 20. |
SAVER |
Whether use and improve SAVER in imputation, default is FALSE. |
MAGIC |
Whether use and improve MAGIC in imputation, default is FALSE. |
UMI |
Whether use full UMI data, default is FALSE. If TRUE, |
hist_raw_counts |
A list contains the histogram table of raw read counts for each cell in |
hist_RUG_counts |
A list contains the histogram table of raw UMI-gene counts for each cell in |
parallel |
If FALSE (default), no parallel computation is used; if TRUE, parallel computation using |
BPPARAM |
An optional parameter object passed internally to |
verbose |
Whether to show specific calculation progress, default is TRUE. |
Imputed counts matrices will be saved in the output directory specified by outputDir
.
Zhun Miao.
estDropoutNum
, for estimating dropout gene number in a cell.
countsSampling
, for downsampling the read counts in a cell.
normalization
, for normalization of single-cell RNA-seq data.
scRecoverTest
, a test dataset for scRecover.
# Load test data for scRecover data(scRecoverTest) # Run scRecover with Kcluster specified scRecover(counts = counts, Kcluster = 2) # Or run scRecover with labels specified # scRecover(counts = counts, labels = labels)
# Load test data for scRecover data(scRecoverTest) # Run scRecover with Kcluster specified scRecover(counts = counts, Kcluster = 2) # Or run scRecover with labels specified # scRecover(counts = counts, labels = labels)
A test dataset containing a single-cell RNA-seq (scRNA-seq) read counts matrix and its cell type information.
data(scRecoverTest)
data(scRecoverTest)
counts. A non-negative integer matrix of scRNA-seq raw read counts, rows are genes and columns are cells.
labels. A vector of integer specifying the cell types in the read counts matrix, corresponding to the columns of counts
.
oneCell. A non-negative vector of scRNA-seq raw read counts of a cell for each gene.
counts. A matrix of raw read counts of scRNA-seq data which has 200 genes (rows) and 150 cells (columns).
labels. A vector of integer specifying the two cell types in counts
. Also could be generated by: labels <- c(rep(1,50), rep(2,100))
.
oneCell. A vector of a cell's raw read counts for 24538 gene.
Petropoulos S, et al. Cell, 2016, 165(4): 1012-1026.
scRecover
, for imputation of single-cell RNA-seq data.
estDropoutNum
, for estimating dropout gene number in a cell.
countsSampling
, for downsampling the read counts in a cell.
normalization
, for normalization of single-cell RNA-seq data.
# Load test data for scRecover data(scRecoverTest)
# Load test data for scRecover data(scRecoverTest)