Contains UMI droplet-based single-cell RNA-seq data.
Usage
CellRanger(
dir,
filtered = TRUE,
organism = NULL,
ensemblRelease = NULL,
genomeBuild = NULL,
gffFile = NULL,
refdataDir = NULL,
samples = NULL,
censorSamples = NULL,
sampleMetadataFile = NULL,
transgeneNames = NULL,
interestingGroups = "sampleName"
)
Arguments
- dir
character(1)
. Directory path to Cell Ranger output.- filtered
logical(1)
. Use filtered (recommended) or raw counts. Note that raw counts still contain only whitelisted cellular barcodes.- organism
character(1)
. Full Latin organism name (e.g."Homo sapiens"
).- ensemblRelease
integer(1)
. Ensembl release version (e.g.100
). We recommend setting this value if possible, for improved reproducibility. When left unset, the latest release available via AnnotationHub/ensembldb is used. Note that the latest version available can vary, depending on the versions of AnnotationHub and ensembldb in use.- genomeBuild
character(1)
. Ensembl genome build assembly name (e.g."GRCh38"
). If setNULL
, defaults to the most recent build available. Note: don't pass in UCSC build IDs (e.g."hg38"
).- gffFile
character(1)
. GFF/GTF (General Feature Format) file. Generally, we recommend using GTF (GFFv2) instead of GFFv3.- refdataDir
character(1)
orNULL
. Directory path to Cell Ranger reference annotation data.- samples
character
. Sample identifiers.- censorSamples
character
. Specify a subset of samples to censor.- sampleMetadataFile
character(1)
. Sample metadata file path. CSV or TSV is preferred, but Excel worksheets are also supported. Check the documentation for conventions and required columns.- transgeneNames
character
. Vector indicating which assay rows denote transgenes (e.g. EGFP, TDTOMATO).- interestingGroups
character
. Groups of interest to use for visualization. Corresponds to factors describing the columns of the object.
Details
Read 10x Genomics Cell Ranger output
for a Chromium data set into a SingleCellExperiment
object.
Currently supports loading of a single genome.
Directory structure for multiple samples
Cell Ranger can vary in its output directory structure, but we're requiring a
single, consistent directory structure for datasets containing multiple
samples that have not been aggregated into a single matrix with aggr
.
Cell Ranger v3 output:
| <dir>/
|-- <sampleName>/
|---- SC_RNA_COUNTER_CS/
|---- outs/
|------ filtered_feature_bc_matrix/
|-------- barcodes.tsv.gz
|-------- features.tsv.gz
|-------- matrix.mtx.gz
|------ filtered_feature_bc_matrix.h5
|------ metrics_summary.csv
|------ molecule_info.h5
|------ possorted_genome_bam.bam
|------ possorted_genome_bam.bam.bai
|------ raw_feature_bc_matrix/
|-------- barcodes.tsv.gz
|-------- features.tsv.gz
|-------- matrix.mtx.gz
|------ raw_feature_bc_matrix.h5
|------ web_summary.html
Cell Ranger v2 output:
| <dir>/
|-- <sampleName>/
|---- SC_RNA_COUNTER_CS/
|---- outs/
|------ filtered_gene_bc_matrices/
|-------- <genomeBuild>/
|---------- barcodes.tsv
|---------- genes.tsv
|---------- matrix.mtx
|------ filtered_gene_bc_matrices_h5.h5
|------ metrics_summary.csv
|------ molecule_info.h5
|------ possorted_genome_bam.bam
|------ possorted_genome_bam.bam.bai
|------ raw_gene_bc_matrices/
|-------- <genomeBuild>/
|---------- barcodes.tsv
|---------- genes.tsv
|---------- matrix.mtx
|------ raw_gene_bc_matrices_h5.h5
Sample metadata
A user-supplied sample metadata file defined by sampleMetadataFile
is
required for multiplexed datasets. Otherwise this can be left NULL
, and
minimal sample data will be used, based on the directory names.
Reference data
We strongly recommend supplying the corresponding reference data required for
Cell Ranger with the refdataDir
argument. It will convert the gene
annotations defined in the GTF file into a GRanges
object, which get
slotted in rowRanges()
. Otherwise, the
function will attempt to use the most current annotations available from
Ensembl, and some gene IDs may not match, due to deprecation in the current
Ensembl release.
Examples
dir <- system.file("extdata", "cellranger_v3", package = "Chromium")
x <- CellRanger(dir)
#> → Importing Chromium single-cell RNA-seq run.
#> ℹ 1 sample detected:
#> • pbmc
#> → Importing counts from matrix.mtx.gz file.
#> → Importing /private/var/folders/l1/8y8sjzmn15v49jgrqglghcfr0000gn/T/RtmpCg1NZ7/temp_libpath15e75769d13a0/Chromium/extdata/cellranger_v3/pbmc/outs/filtered_feature_bc_matrix/matrix.mtx.gz using Matrix::`readMM()`.
#> → Importing /private/var/folders/l1/8y8sjzmn15v49jgrqglghcfr0000gn/T/RtmpCg1NZ7/temp_libpath15e75769d13a0/Chromium/extdata/cellranger_v3/pbmc/outs/filtered_feature_bc_matrix/features.tsv.gz using base::`read.table()`.
#> → Importing /private/var/folders/l1/8y8sjzmn15v49jgrqglghcfr0000gn/T/RtmpCg1NZ7/temp_libpath15e75769d13a0/Chromium/extdata/cellranger_v3/pbmc/outs/filtered_feature_bc_matrix/barcodes.tsv.gz using base::`read.table()`.
#> ! Slotting empty ranges into `rowRanges()`.
#> ℹ Filtered zero count rows and columns:
#> - 46 / 100 rows (46%)
#> - 92 / 100 columns (92%)
#> → Calculating 92 sample metrics.
#> ! Calculating metrics without biotype information.
#> `rowData()` required to calculate: "nCoding", "nMito", "mitoRatio".
#> ℹ Prefilter: 55 / 92 samples (60%).
#> ✔ Chromium single-cell RNA-seq run imported successfully.
print(x)
#> class: CellRanger
#> dim: 46 55
#> metadata(21): allSamples call ... umiType wd
#> assays(1): counts
#> rownames(46): ENSG00000008128 ENSG00000008130 ... ENSG00000272106
#> ENSG00000272512
#> rowData names(0):
#> colnames(55): AAACCCAAGGCCTAGA AAACCCATCGTGCATA ... AACCTGAGTCGAGCTC
#> AACCTGAGTCTGCAAT
#> colData names(8): sampleId sampleName ... log10FeaturesPerCount
#> mitoRatio
#> reducedDimNames(0):
#> mainExpName: NULL
#> altExpNames(0):