QuickCB2 | R Documentation |
All-in-one function for scCB2
. Take 10x output raw data as input
and return either a matrix of real cells identified by CB2 or
a Seurat object containing this matrix, which can be incorporated
with downstream analysis using Seurat pipeline.
QuickCB2(
dir = NULL,
h5file = NULL,
FDR_threshold = 0.01,
MTfilter = 1,
MTgene = NULL,
AsSeurat = FALSE,
Ncores = 2,
...
)
dir |
The directory of 10x output data. For Cell Ranger version <3, directory should include three files: barcodes.tsv, genes.tsv, matrix.mtx. For Cell Ranger version >=3, directory should include three files: barcodes.tsv.gz, features.tsv.gz, matrix.mtx.gz. |
h5file |
The path of 10x output HDF5 file (ended with .h5). |
FDR_threshold |
Numeric between 0 and 1. Default: 0.01. The False Discovery Rate (FDR) to be controlled for multiple testing. |
MTfilter |
Numeric value between 0 and 1. Default: |
MTgene |
Character vector. User may specify customized mitochondrial gene IDs to perform the filtering. This should correspond to a subset of row names in raw data. |
AsSeurat |
Logical. Default: |
Ncores |
Positive integer. Default: |
... |
Additional arguments to be passed to |
QuickCB2
is a quick function to apply CB2 on 10x Cell Ranger
raw data by combining Read10xRaw
, Read10xRawH5
,
CB2FindCell
and GetCellMat
into one simple function
under default parameters.
Either a sparse matrix of real cells identified by CB2 or a Seurat object containing real cell matrix.
# simulate 10x output files
data(mbrainSub)
mbrainSub <- mbrainSub[,1:10000]
data_dir <- file.path(tempdir(),"CB2example")
dir.create(data_dir)
gene_name <- rownames(mbrainSub)
# For simplicity, use gene names to generate gene IDs to fit the format.
gene_id <- paste0("ENSG_fake_",gene_name)
barcode_id <- colnames(mbrainSub)
Matrix::writeMM(mbrainSub,file = file.path(data_dir,"matrix.mtx"))
write.table(barcode_id,file = file.path(data_dir,"barcodes.tsv"),
sep = "\t", quote = FALSE, col.names = FALSE, row.names = FALSE)
write.table(cbind(gene_id,gene_name),file = file.path(data_dir,"genes.tsv"),
sep = "\t", quote = FALSE, col.names = FALSE, row.names = FALSE)
# Run QuickCB2 on 10x raw data and get cell matrix.
# Control FDR at 1%. Use 2-core parallel computation.
RealCell <- QuickCB2(dir = data_dir,
FDR_threshold = 0.01,
Ncores = 2)
str(RealCell)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.