%\VignetteIndexEntry{Counting reads with summarizeOverlaps} %\VignetteDepends{GenomicAlignments, DESeq2, edgeR, rtracklayer} %\VignetteKeywords{sequence, sequencing, alignments} %\VignettePackage{GenomicAlignments} \documentclass{article} <>= BiocStyle::latex() @ \title{Counting reads with \Rfunction{summarizeOverlaps}} \author{Valerie Obenchain} \date{Edited: February 2024; Compiled: \today} \begin{document} \maketitle \tableofcontents <>= options(width=72) options("showHeadLines" = 3) options("showTailLines" = 3) @ \section{Introduction} This vignette illustrates how reads mapped to a genome can be counted with \Rfunction{summarizeOverlaps}. Different "modes" of counting are provided to resolve reads that overlap multiple features. The built-in count modes are fashioned after the "Union", "IntersectionStrict", and "IntersectionNotEmpty" methods found in the HTSeq package by Simon Anders (see references). \section{A First Example} In this example reads are counted from a list of BAM files and returned in a \Robject{matrix} for use in further analysis such as those offered in \Biocpkg{DESeq2} and \Biocpkg{edgeR}. <>= library(GenomicAlignments) library(DESeq2) library(edgeR) fls <- list.files(system.file("extdata", package="GenomicAlignments"), recursive=TRUE, pattern="*bam$", full=TRUE) features <- GRanges( seqnames = c(rep("chr2L", 4), rep("chr2R", 5), rep("chr3L", 2)), ranges = IRanges(c(1000, 3000, 4000, 7000, 2000, 3000, 3600, 4000, 7500, 5000, 5400), width=c(rep(500, 3), 600, 900, 500, 300, 900, 300, 500, 500)), "-", group_id=c(rep("A", 4), rep("B", 5), rep("C", 2))) olap <- summarizeOverlaps(features, fls) deseq <- DESeqDataSet(olap, design= ~ 1) edger <- DGEList(assay(olap), group=rownames(colData(olap))) @ %% By default, the \Rfunction{summarizeOverlaps} function iterates through files in `chunks' and with files processed in parallel. For finer-grain control over memory consumption, use the \Rfunction{BamFileList} function and specify the \Rcode{yieldSize} argument (e.g., \Rcode{yieldSize=1000000}) to determine the size of each `chunk' (smaller chunks consume less memory, but are a little less efficient to process). For controlling the number of processors in use, use \Rfunction{BiocParallel::register} to use an appropriate back-end, e.g., in linux or Mac to process on 6 cores of a single machine use \Rcode{register(MulticoreParam(workers=6))}; see the \Biocpkg{BiocParallel} vignette for further details. \section{Counting Modes} The modes of "Union", "IntersectionStrict" and "IntersectionNotEmpty" provide different approaches to resolving reads that overlap multiple features. Figure~\ref{fig-summarizeOverlaps-modes} illustrates how both simple and gapped reads are handled by the modes. Note that a read is counted a maximum of once; there is no double counting. For additional detail on the counting modes see the \Rfunction{summarizeOverlaps} man page. \begin{figure}[!h] \begin{center} \includegraphics{summarizeOverlaps-modes.pdf} \caption{Counting Modes} \label{fig-summarizeOverlaps-modes} \end{center} \end{figure} \newpage \section{Counting Features} Features can be exons, transcripts, genes or any region of interest. The number of ranges that define a single feature is specified in the \Rcode{features} argument. When annotation regions of interest are defined by a single range a \Rclass{GRanges} should be used as the \Rcode{features} argument. With a \Rclass{GRanges} it is assumed that each row (i.e., each range) represents a distinct feature. If \Rcode{features} was a \Rclass{GRanges} of exons, the result would be counts per exon. When the region of interest is defined by one or more ranges the \Rcode{features} argument should be a \Rclass{GRangesList}. In practice this could be a list of exons by gene or transcripts by gene or other similar relationships. The count result will be the same length as the \Rclass{GRangesList}. For a list of exons by genes, the result would be counts per gene. The combination of defining the features as either\Rclass{GRanges} or \Rclass{GRangesList} and choosing a counting mode controls how \Rfunction{summarizeOverlaps} assigns hits. Regardless of the mode chosen, each read is assigned to at most a single feature. These options are intended to provide flexibility in defining different biological problems. This next example demonstrates how the same read can be counted differently depending on how the \Rcode{features} argument is specified. We use a single read that overlaps two ranges, gr1 and gr2. <>= rd <- GAlignments(names="a", seqnames="chr1", pos=100, cigar="300M", strand="+") gr1 <- GRanges("chr1", IRanges(start=50, width=150), strand="+") gr2 <- GRanges("chr1", IRanges(start=350, width=150), strand="+") @ \noindent When provided as a \Rclass{GRanges} both gr1 and gr2 are considered distinct features. In this case none of the modes count the read as a hit. Mode \Rcode{Union} discards the read becasue more than 1 feature is overlapped. \Rcode{IntersectionStrict} requires the read to fall completely within a feature which is not the case for either gr1 or gr2. \Rcode{IntersetctionNotEmpty} requires the read to overlap a single unique disjoint region of the \Rcode{features}. In this case gr1 and gr2 do not overlap so each range is considered a unique disjoint region. However, the read overlaps both gr1 and gr2 so a decision cannot be made and the read is discarded. <>= gr <- c(gr1, gr2) data.frame(union = assay(summarizeOverlaps(gr, rd)), intStrict = assay(summarizeOverlaps(gr, rd, mode="IntersectionStrict")), intNotEmpty = assay(summarizeOverlaps(gr, rd, mode="IntersectionNotEmpty"))) @ \noindent Next we count with \Rcode{features} as a \Rclass{GRangesList}; this is list of length 1 with 2 elements. Modes \Rcode{Union} and \Rcode{IntersectionNotEmpty} both count the read for the single feature. <>= grl <- GRangesList(c(gr1, gr2)) data.frame(union = assay(summarizeOverlaps(grl, rd)), intStrict = assay(summarizeOverlaps(grl, rd, mode="IntersectionStrict")), intNotEmpty = assay(summarizeOverlaps(grl, rd, mode="IntersectionNotEmpty"))) @ In this more complicated example we have 7 reads, 5 are simple and 2 have gaps in the CIGAR. There are 12 ranges that will serve as the \Robject{features}. <>= group_id <- c("A", "B", "C", "C", "D", "D", "E", "F", "G", "G", "H", "H") features <- GRanges( seqnames = c("chr1", "chr2", "chr1", "chr1", "chr2", "chr2", "chr1", "chr1", "chr2", "chr2", "chr1", "chr1"), strand = Rle(strand("+"), length(group_id)), ranges = IRanges( start=c(1000, 2000, 3000, 3600, 7000, 7500, 4000, 4000, 3000, 3350, 5000, 5400), width=c(500, 900, 500, 300, 600, 300, 500, 900, 150, 200, 500, 500)), DataFrame(group_id) ) reads <- GAlignments( names = c("a","b","c","d","e","f","g"), seqnames = Rle(c(rep(c("chr1", "chr2"), 3), "chr1")), pos = c(1400, 2700, 3400, 7100, 4000, 3100, 5200), cigar = c("500M", "100M", "300M", "500M", "300M", "50M200N50M", "50M150N50M"), strand = Rle(strand("+"), 7)) @ \noindent Using a \Rclass{GRanges} as the \Rcode{features} all 12 ranges are considered to be different features and counts are produced for each row, <>= data.frame(union = assay(summarizeOverlaps(features, reads)), intStrict = assay(summarizeOverlaps(features, reads, mode="IntersectionStrict")), intNotEmpty = assay(summarizeOverlaps(features, reads, mode="IntersectionNotEmpty"))) @ \noindent When the data are split by group to create a \Rclass{GRangesList} the highest list-levels are treated as different features and the multiple list elements are considered part of the same features. Counts are returned for each group. <>= lst <- split(features, mcols(features)[["group_id"]]) length(lst) @ <>= data.frame(union = assay(summarizeOverlaps(lst, reads)), intStrict = assay(summarizeOverlaps(lst, reads, mode="IntersectionStrict")), intNotEmpty = assay(summarizeOverlaps(lst, reads, mode="IntersectionNotEmpty"))) @ If desired, users can supply their own counting function as the \Rcode{mode} argument and take advantage of the infrastructure for counting over multiple BAM files and parsing the results into a \Rclass{RangedSummarizedExperiment} object. See \Rcode{?'BamViews-class'} or \Rcode{?'BamFile-class'} in the \Biocpkg{Rsamtools} package. \section{\Rcode{pasilla} Data} In this excercise we count the \Biocpkg{pasilla} data by gene and by transcript then create a \Rclass{DESeqDataSet}. This object can be used in differential expression methods offered in the \Biocpkg{DESeq2} package. \subsection{source files} Files are available through NCBI Gene Expression Omnibus (GEO), accession number GSE18508. \url{http://www.ncbi.nlm.nih.gov/projects/geo/query/acc.cgi?acc=GSE18508}. SAM files can be converted to BAM with the \Rfunction{asBam} function in the \Biocpkg{Rsamtools} package. Of the seven files available, 3 are single-reads and 4 are paired-end. Smaller versions of untreated1 (single-end) and untreated2 (paired-end) have been made available in the \Biocpkg{pasillaBamSubset} package. This subset includes chromosome 4 only. \Rfunction{summarizeOverlaps} is capable of counting paired-end reads in both a \Rcode{BamFile}-method (set argument \Rcode{singleEnd=TRUE}) or a \Rcode{GAlignmentPairs}-method. For this example, we use the 3 single-end read files, \begin{itemize} \item treated1.bam \item untreated1.bam \item untreated2.bam \end{itemize} Annotations are retrieved as a GTF file from the ENSEMBL web site. We download the file our local disk, then use \Biocpkg{Rtracklayer}'s \Rfunction{import} function to parse the file to a \Rclass{GRanges} instance. <>= library(rtracklayer) fl <- paste0("ftp://ftp.ensembl.org/pub/release-62/", "gtf/drosophila_melanogaster/", "Drosophila_melanogaster.BDGP5.25.62.gtf.gz") gffFile <- file.path(tempdir(), basename(fl)) download.file(fl, gffFile) gff0 <- import(gffFile) @ Subset on the protein-coding, exon regions of chromosome 4 and split by gene id. <>= idx <- mcols(gff0)$source == "protein_coding" & mcols(gff0)$type == "exon" & seqnames(gff0) == "4" gff <- gff0[idx] ## adjust seqnames to match Bam files seqlevels(gff) <- paste("chr", seqlevels(gff), sep="") chr4genes <- split(gff, mcols(gff)$gene_id) @ \subsection{counting} The \Rcode{param} argument can be used to subset the reads in the bam file on characteristics such as position, unmapped or paired-end reads. Quality scores or the "NH" tag, which identifies reads with multiple mappings, can be included as metadata columns for further subsetting. See \Rcode{?ScanBamParam} for details about specifying the \Rcode{param} argument. <>= param <- ScanBamParam( what='qual', which=GRanges("chr4", IRanges(1, 1e6)), flag=scanBamFlag(isUnmappedQuery=FALSE, isPaired=NA), tag="NH") @ We use \Rfunction{summarizeOverlaps} to count with the default mode of "Union". If a \Rcode{param} argument is not included all reads from the BAM file are counted. <>= fls <- c("treated1.bam", "untreated1.bam", "untreated2.bam") path <- "pathToBAMFiles" bamlst <- BamFileList(fls) genehits <- summarizeOverlaps(chr4genes, bamlst, mode="Union") @ \noindent A \Rcode{CountDataSet} is constructed from the counts and experiment data in \Rclass{pasilla}. <>= expdata <- MIAME( name="pasilla knockdown", lab="Genetics and Developmental Biology, University of Connecticut Health Center", contact="Dr. Brenton Graveley", title="modENCODE Drosophila pasilla RNA Binding Protein RNAi knockdown RNA-Seq Studies", pubMedIds="20921232", url="http://www.ncbi.nlm.nih.gov/projects/geo/query/acc.cgi?acc=GSE18508", abstract="RNA-seq of 3 biological replicates of from the Drosophila melanogaster S2-DRSC cells that have been RNAi depleted of mRNAs encoding pasilla, a mRNA binding protein and 4 biological replicates of the the untreated cell line.") design <- data.frame( condition=c("treated", "untreated", "untreated"), replicate=c(1,1,2), type=rep("single-read", 3), countfiles=path(colData(genehits)[,1]), stringsAsFactors=TRUE) geneCDS <- DESeqDataSet(genehits, design=design, metadata=list(expdata=expdata)) @ If the primary interest is to count by transcript instead of by gene, the annotation file can be split on transcript id. <>= chr4tx <- split(gff, mcols(gff)$transcript_id) txhits <- summarizeOverlaps(chr4tx, bamlst) txCDS <- DESeqDataSet(txhits, design=design, metadata=list(expdata=expdata)) @ \section{References} \url{http://www-huber.embl.de/users/anders/HTSeq/doc/overview.html} \noindent\url{http://www-huber.embl.de/users/anders/HTSeq/doc/count.html} \end{document}