View source: R/importCellRanger.R
importCellRanger | R Documentation |
Read the filtered barcodes, features, and matrices for all samples from (preferably a single run of) Cell Ranger output. Import and combine them as one big SingleCellExperiment object.
importCellRanger(
cellRangerDirs = NULL,
sampleDirs = NULL,
sampleNames = NULL,
cellRangerOuts = NULL,
dataType = c("filtered", "raw"),
matrixFileNames = "matrix.mtx.gz",
featuresFileNames = "features.tsv.gz",
barcodesFileNames = "barcodes.tsv.gz",
gzipped = "auto",
class = c("Matrix", "matrix"),
delayedArray = FALSE,
rowNamesDedup = TRUE
)
importCellRangerV2(
cellRangerDirs = NULL,
sampleDirs = NULL,
sampleNames = NULL,
dataTypeV2 = c("filtered", "raw"),
class = c("Matrix", "matrix"),
delayedArray = FALSE,
reference = NULL,
cellRangerOutsV2 = NULL,
rowNamesDedup = TRUE
)
importCellRangerV3(
cellRangerDirs = NULL,
sampleDirs = NULL,
sampleNames = NULL,
dataType = c("filtered", "raw"),
class = c("Matrix", "matrix"),
delayedArray = FALSE,
rowNamesDedup = TRUE
)
cellRangerDirs |
The root directories where Cell Ranger was run. These
folders should contain sample specific folders. Default |
sampleDirs |
Default
The cells in the final SCE object will be ordered in the same order of
|
sampleNames |
A vector of user-defined sample names for the samples
to be
imported. Must have the same length as |
cellRangerOuts |
Character vector. The intermediate
paths to filtered or raw cell barcode, feature, and matrix files
for each sample. Supercedes |
dataType |
Character. The type of data to import. Can be one of
"filtered" (which is equivalent to
|
matrixFileNames |
Character vector. Filenames for the Market Exchange
Format (MEX) sparse matrix files (matrix.mtx or matrix.mtx.gz files).
Must have length 1 or the same
length as |
featuresFileNames |
Character vector. Filenames for the feature
annotation files. They are usually named features.tsv.gz or
genes.tsv. Must have length 1 or the same
length as |
barcodesFileNames |
Character vector. Filename for the cell barcode
list files. They are usually named barcodes.tsv.gz or
barcodes.tsv. Must have length 1 or the same
length as |
gzipped |
|
class |
Character. The class of the expression matrix stored in the SCE
object. Can be one of "Matrix" (as returned by
readMM function), or "matrix" (as returned by
matrix function). Default |
delayedArray |
Boolean. Whether to read the expression matrix as
DelayedArray object or not. Default |
rowNamesDedup |
Boolean. Whether to deduplicate rownames. Default
|
dataTypeV2 |
Character. The type of output to import for
Cellranger version below 3.0.0. Whether to import the filtered or the
raw data. Can be one of 'filtered' or 'raw'. Default 'filtered'. When
|
reference |
Character vector. The reference genome names.
Default |
cellRangerOutsV2 |
Character vector. The intermediate paths
to filtered or raw cell barcode, feature, and matrix files for each
sample for Cellranger version below 3.0.0. If |
importCellRangerV2
imports output from Cell Ranger V2.
importCellRangerV2Sample
imports output from one sample from Cell
Ranger V2.
importCellRangerV3
imports output from Cell Ranger V3.
importCellRangerV3
imports output from one sample from Cell Ranger
V3.
Some implicit
assumptions which match the output structure of Cell Ranger V2 & V3
are made in these 4 functions including cellRangerOuts
,
matrixFileName
, featuresFileName
, barcodesFileName
,
and gzipped
.
Alternatively, user can call importCellRanger
to explicitly
specify these arguments.
A SingleCellExperiment
object containing the combined count
matrix, the feature annotations, and the cell annotation.
# Example #1
# The following filtered feature, cell, and matrix files were downloaded from
# https://support.10xgenomics.com/single-cell-gene-expression/datasets/
# 3.0.0/hgmm_1k_v3
# The top 10 hg19 & mm10 genes are included in this example.
# Only the first 20 cells are included.
sce <- importCellRanger(
cellRangerDirs = system.file("extdata/", package = "singleCellTK"),
sampleDirs = "hgmm_1k_v3_20x20",
sampleNames = "hgmm1kv3",
dataType = "filtered")
# The following filtered feature, cell, and matrix files were downloaded from
# https://support.10xgenomics.com/single-cell-gene-expression/datasets/
# 2.1.0/pbmc4k
# Top 20 genes are kept. 20 cell barcodes are extracted.
sce <- importCellRangerV2(
cellRangerDirs = system.file("extdata/", package = "singleCellTK"),
sampleDirs = "pbmc_4k_v2_20x20",
sampleNames = "pbmc4k_20",
reference = 'GRCh38',
dataTypeV2 = "filtered")
sce <- importCellRangerV3(
cellRangerDirs = system.file("extdata/", package = "singleCellTK"),
sampleDirs = "hgmm_1k_v3_20x20",
sampleNames = "hgmm1kv3",
dataType = "filtered")
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.