tenxBamqc | R Documentation |
Read BAM file generated by Cell Ranger pipeline and output QC metrics including number of aligned reads and reads aligned to an gene.
tenxBamqc(
bam,
experiment,
filter,
validCb = NA,
tags = c("NH", "GX", "CB", "MM"),
yieldSize = 1e+06,
outDir = "./",
cores = max(1, parallelly::availableCores() - 2)
)
bam |
Paths to input BAM files generated by Cell Ranger pipeline. These files are usually named "possorted_genome_bam.bam" in the "outs" folder of the top-level project output folders, respectively. |
experiment |
A character vector of experiment names. Represents the
group label for each BAM file, e.g. "patient1, patient2, ...". The
length of |
filter |
Paths to the filtered barcode files. Should be in the same length and order of the input BAM files. These files are named "barcodes.tsv" located at outs/filtered_gene_bc_matrices/<reference_genome>/barcodes.tsv. |
validCb |
Path to the cell barcode whitelist file. By default uses the
file "737K-august-2016.txt" which is compatible with the v2 chemistry
protocol. The file can be inspected by calling
|
tags |
BAM tags used for collecting QC metrics. Contains non-standard tags locally-defined by Cell Ranger pipeline. Should not be changed in most cases. |
yieldSize |
The number of records (alignments) to yield when drawing
successive subsets from a BAM file, providing the number of successive
records to be returned on each yield. This parameter is passed to the
|
outDir |
Output directory. The location to write resulting QC table. |
cores |
Number of cores used for parallelization. Default is
|
A SingleCellExperiment object. The
colData
contains the number of aligned reads
(reads_mapped_to_genome) and reads aligned
to genes (reads_mapped_to_genes).
# first 5000 records in the bam file downloaded from here:
# http://sra-download.ncbi.nlm.nih.gov/srapub_files/
# SRR5167880_E18_20160930_Neurons_Sample_01.bam
# see details here:
# https://trace.ncbi.nlm.nih.gov/Traces/sra/?study=SRP096558
# and here:
# https://www.ncbi.nlm.nih.gov/geo/query/acc.cgi?acc=GSE93421
bamfile10x <- system.file("extdata",
"SRR5167880_E18_20160930_Neurons_Sample_01_5000.bam",
package = "scruff")
# library(TENxBrainData)
# library(data.table)
# tenx <- TENxBrainData()
# # get filtered barcodes for sample 01
# filteredBcIndex <- tstrsplit(colData(tenx)[, "Barcode"], "-")[[2]] == 1
# filteredBc <- colData(tenx)[filteredBcIndex, ][["Barcode"]]
filteredBc <- system.file("extdata",
"SRR5167880_E18_20160930_Neurons_Sample_01_filtered_barcode.tsv",
package = "scruff")
# QC results are saved to current working directory
sce <- tenxBamqc(bam = bamfile10x,
experiment = "Neurons_Sample_01",
filter = filteredBc)
sce
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.