%\VignetteIndexEntry{2. GenomicRanges HOWTOs} %\VignetteDepends{GenomicRanges, Rsamtools, GenomicAlignments, pasillaBamSubset, TxDb.Dmelanogaster.UCSC.dm3.ensGene, TxDb.Athaliana.BioMart.plantsmart22, AnnotationHub, DESeq2, edgeR, TxDb.Hsapiens.UCSC.hg38.knownGene, GenomicFeatures, txdbmaker, Biostrings, BSgenome.Hsapiens.UCSC.hg38, KEGGREST, KEGGgraph, BSgenome.Scerevisiae.UCSC.sacCer2} %\VignetteKeywords{sequence, sequencing, alignments} %\VignettePackage{GenomicRanges} \documentclass{article} \usepackage[authoryear,round]{natbib} <>= BiocStyle::latex(use.unsrturl=FALSE) @ \title{\Biocpkg{GenomicRanges} HOWTOs} \author{Bioconductor Team} \date{Edited: June 2024; Compiled: \today} \begin{document} \maketitle \tableofcontents <>= options(width=72) options(showHeadLines=3) options(showTailLines=3) .precomputed_results_path <- "precomputed_results" @ %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% \section{Introduction} \subsection{Purpose of this document} This document is a collection of {\it HOWTOs}. Each {\it HOWTO} is a short section that demonstrates how to use the containers and operations implemented in the \Biocpkg{GenomicRanges} and related packages (\Biocpkg{IRanges}, \Biocpkg{Biostrings}, \Biocpkg{Rsamtools}, \Biocpkg{GenomicAlignments}, \Biocpkg{BSgenome}, and \Biocpkg{GenomicFeatures}) to perform a task typically found in the context of a high throughput sequence analysis. Unless stated otherwise, the {\it HOWTOs} are self contained, independent of each other, and can be studied and reproduced in any order. \subsection{Prerequisites and additional recommended reading} We assume the reader has some previous experience with \R{} and with basic manipulation of \Rcode{GRanges}, \Rcode{GRangesList}, \Rcode{Rle}, \Rcode{RleList}, and \Rcode{DataFrame} objects. See the ``An Introduction to Genomic Ranges Classes'' vignette located in the \Biocpkg{GenomicRanges} package (in the same folder as this document) for an introduction to these containers. Additional recommended readings after this document are the ``Software for Computing and Annotating Genomic Ranges'' paper[\citet{Lawrence2013ranges}] and the ``Counting reads with \Rfunction{summarizeOverlaps}'' vignette located in the \Biocpkg{GenomicAlignments} package. To display the list of vignettes available in the \Biocpkg{GenomicRanges} package, use \Rcode{browseVignettes("GenomicRanges")}. \subsection{Input data and terminology used across the HOWTOs} In order to avoid repetition, input data, concepts and terms used in more than one {\it HOWTO} are described here: \begin{itemize} \item {\bf The \Biocpkg{pasillaBamSubset} data package}: contains both a BAM file with single-end reads (untreated1\_chr4) and a BAM file with paired-end reads (untreated3\_chr4). Each file is a subset of chr4 from the "Pasilla" experiment. <>= library(pasillaBamSubset) untreated1_chr4() untreated3_chr4() @ See \Rcode{?pasillaBamSubset} for more information. <>= ?pasillaBamSubset @ \item {\bf Gene models and \Rclass{TxDb} objects}: A \textit{gene model} is essentially a set of annotations that describes the genomic locations of the known genes, transcripts, exons, and CDS, for a given organism. In \Bioconductor{} it is typically represented as a \Rclass{TxDb} object but also sometimes as a \Rclass{GRanges} or \Rclass{GRangesList} object. The \Biocpkg{txdbmaker} and \Biocpkg{GenomicFeatures} packages contain tools for making and manipulating \Rclass{TxDb} objects. \end{itemize} %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% \section{HOWTOs} %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% \subsection{How to read single-end reads from a BAM file} As sample data we use the \Biocpkg{pasillaBamSubset} data package described in the introduction. <>= library(pasillaBamSubset) un1 <- untreated1_chr4() # single-end reads @ Several functions are available for reading BAM files into \R{}: \begin{verbatim} readGAlignments() readGAlignmentPairs() readGAlignmentsList() scanBam() \end{verbatim} \Rfunction{scanBam} is a low-level function that returns a list of lists and is not discussed further here. See \Rcode{?scanBam} in the \Biocpkg{Rsamtools} package for more information. Single-end reads can be loaded with the \Rfunction{readGAlignments} function from the \Biocpkg{GenomicAlignments} package. <>= library(GenomicAlignments) gal <- readGAlignments(un1) @ Data subsets can be specified by genomic position, field names, or flag criteria in the \Rcode{ScanBamParam}. Here we input records that overlap position 1 to 5000 on the negative strand with \Rcode{flag} and \Rcode{cigar} as metadata columns. <>= what <- c("flag", "cigar") which <- GRanges("chr4", IRanges(1, 5000)) flag <- scanBamFlag(isMinusStrand = TRUE) param <- ScanBamParam(which=which, what=what, flag=flag) neg <- readGAlignments(un1, param=param) neg @ Another approach to subsetting the data is to use \Rfunction{filterBam}. This function creates a new BAM file of records passing user-defined criteria. See \Rcode{?filterBam} in the \Biocpkg{Rsamtools} package for more information. %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% \subsection{How to read paired-end reads from a BAM file} As sample data we use the \Biocpkg{pasillaBamSubset} data package described in the introduction. <>= library(pasillaBamSubset) un3 <- untreated3_chr4() # paired-end reads @ Paired-end reads can be loaded with the \Rfunction{readGAlignmentPairs} or \Rfunction{readGAlignmentsList} function from the \Biocpkg{GenomicAlignments} package. These functions use the same mate paring algorithm but output different objects. Let's start with \Rfunction{readGAlignmentPairs}: <>= un3 <- untreated3_chr4() gapairs <- readGAlignmentPairs(un3) @ The \Robject{GAlignmentPairs} class holds only pairs; reads with no mate or with ambiguous pairing are discarded. Each list element holds exactly 2 records (a mated pair). Records can be accessed as the \Rcode{first} and\Rcode{last} segments in a template or as \Rcode{left} and \Rcode{right} alignments. See \Rcode{?GAlignmentPairs} in the \Biocpkg{GenomicAlignments} package for more information. <>= gapairs @ For \Rcode{readGAlignmentsList}, mate pairing is performed when \Rcode{asMates} is set to \Rcode{TRUE} on the \Rcode{BamFile} object, otherwise records are treated as single-end. <>= galist <- readGAlignmentsList(BamFile(un3, asMates=TRUE)) @ \Robject{GAlignmentsList} is a more general `list-like' structure that holds mate pairs as well as non-mates (i.e., singletons, records with unmapped mates etc.) A \Rcode{mates\_status} metadata column (accessed with \Rfunction{mcols}) indicates which records were paired. <>= galist @ Non-mated reads are returned as groups by QNAME and contain any number of records. Here the non-mate groups range in size from 1 to 9. <>= non_mates <- galist[unlist(mcols(galist)$mate_status) == "unmated"] table(elementNROWS(non_mates)) @ %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% \subsection{How to read and process a big BAM file by chunks in order to reduce memory usage} A large BAM file can be iterated through in chunks by setting a \Rcode{yieldSize} on the \Rclass{BamFile} object. As sample data we use the \Biocpkg{pasillaBamSubset} data package described in the introduction. <>= library(pasillaBamSubset) un1 <- untreated1_chr4() bf <- BamFile(un1, yieldSize=100000) @ Iteration through a BAM file requires that the file be opened, repeatedly queried inside a loop, then closed. Repeated calls to \Rfunction{readGAlignments} without opening the file first result in the same 100000 records returned each time. <>= open(bf) cvg <- NULL repeat { chunk <- readGAlignments(bf) if (length(chunk) == 0L) break chunk_cvg <- coverage(chunk) if (is.null(cvg)) { cvg <- chunk_cvg } else { cvg <- cvg + chunk_cvg } } close(bf) cvg @ %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% \subsection{How to compute read coverage} The ``read coverage'' is the number of reads that cover a given genomic position. Computing the read coverage generally consists in computing the coverage at each position in the genome. This can be done with the \Rcode{coverage()} function. As sample data we use the \Biocpkg{pasillaBamSubset} data package described in the introduction. <>= library(pasillaBamSubset) un1 <- untreated1_chr4() # single-end reads library(GenomicAlignments) reads1 <- readGAlignments(un1) cvg1 <- coverage(reads1) cvg1 @ Coverage on chr4: <>= cvg1$chr4 @ Average and max coverage: <>= mean(cvg1$chr4) max(cvg1$chr4) @ Note that \Rcode{coverage()} is a generic function with methods for different types of objects. See \Rcode{?coverage} for more information. %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% \subsection{How to find peaks in read coverage} ChIP-Seq analysis usually involves finding peaks in read coverage. This process is sometimes called ``peak calling'' or ``peak detection''. Here we're only showing a naive way to find peaks in the object returned by the \Rcode{coverage()} function. \Bioconductor{} packages \Biocpkg{BayesPeak}, \Biocpkg{bumphunter}, \Biocpkg{Starr}, \Biocpkg{CexoR}, \Biocpkg{exomePeak}, \Biocpkg{RIPSeeker}, and others, provide sophisticated peak calling tools for ChIP-Seq, RIP-Seq, and other kind of high throughput sequencing data. Let's assume \Rcode{cvg1} is the object returned by \Rcode{coverage()} (see previous {\it HOWTO} for how to compute it). We can use the \Rcode{slice()} function to find the genomic regions where the coverage is greater or equal to a given threshold. <>= chr4_peaks <- slice(cvg1$chr4, lower=500) chr4_peaks length(chr4_peaks) # nb of peaks @ The weight of a given peak can be defined as the number of aligned nucleotides that belong to the peak (a.k.a. the area under the peak in mathematics). It can be obtained with \Rcode{sum()}: <>= sum(chr4_peaks) @ %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% \subsection{How to retrieve a gene model from the UCSC genome browser} See introduction for a quick description of what \textit{gene models} and \Rclass{TxDb} objects are. We can use the \Rcode{make\-Transcript\-Db\-From\-UCSC()} function from the \Biocpkg{txdbmaker} package to import a UCSC genome browser track as a \Rclass{TxDb} object. <>= library(txdbmaker) ### Internet connection required! Can take several minutes... txdb <- makeTxDbFromUCSC(genome="sacCer2", tablename="ensGene") @ See \Rcode{?makeTxDbFromUCSC} in the \Biocpkg{txdbmaker} package for more information. Note that some of the most frequently used gene models are available as TxDb packages. A TxDb package consists of a pre-made \Rclass{TxDb} object wrapped into an annotation data package. Go to \url{http://bioconductor.org/packages/release/BiocViews.html#\_\_\_TxDb} to browse the list of available TxDb packages. <>= library(TxDb.Hsapiens.UCSC.hg38.knownGene) txdb <- TxDb.Hsapiens.UCSC.hg38.knownGene txdb @ Extract the transcript coordinates from this gene model: <>= transcripts(txdb) @ %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% \subsection{How to retrieve a gene model from Ensembl} See introduction for a quick description of what \textit{gene models} and \Rclass{TxDb} objects are. We can use the \Rcode{make\-Transcript\-Db\-From\-Biomart()} function from the \Biocpkg{txdbmaker} package to retrieve a gene model from the Ensembl Mart. <>= library(txdbmaker) ### Internet connection required! Can take several minutes... txdb <- makeTxDbFromBiomart(biomart="ensembl", dataset="hsapiens_gene_ensembl") @ See \Rcode{?makeTxDbFromBiomart} in the \Biocpkg{txdbmaker} package for more information. Note that some of the most frequently used gene models are available as TxDb packages. A TxDb package consists of a pre-made \Rclass{TxDb} object wrapped into an annotation data package. Go to \url{http://bioconductor.org/packages/release/BiocViews.html#\_\_\_TxDb} to browse the list of available TxDb packages. <>= library(TxDb.Athaliana.BioMart.plantsmart22) txdb <- TxDb.Athaliana.BioMart.plantsmart22 txdb @ Extract the exon coordinates from this gene model: <>= exons(txdb) @ %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% \subsection{How to load a gene model from a GFF or GTF file} See introduction for a quick description of what \textit{gene models} and \Rclass{TxDb} objects are. We can use the \Rcode{make\-Transcript\-Db\-From\-GFF()} function from the \Biocpkg{txdbmaker} package to import a GFF or GTF file as a \Rclass{TxDb} object. <>= library(txdbmaker) gff_file <- system.file("extdata", "GFF3_files", "a.gff3", package="txdbmaker") txdb <- makeTxDbFromGFF(gff_file, format="gff3") txdb @ See \Rcode{?makeTxDbFromGFF} in the \Biocpkg{txdbmaker} package for more information. Extract the exon coordinates grouped by gene from this gene model: <>= exonsBy(txdb, by="gene") @ %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% \subsection{How to retrieve a gene model from \Biocpkg{AnnotationHub}} When a gene model is not available as a \Rclass{GRanges} or \Rclass{GRangesList} object or as a \Bioconductor{} data package, it may be available on \Biocpkg{AnnotationHub}. In this {\it HOWTO}, will look for a gene model for Drosophila melanogaster on \Biocpkg{AnnotationHub}. Create a `hub' and then filter on Drosophila melanogaster: <>= library(AnnotationHub) ### Internet connection required! hub <- AnnotationHub() hub <- subset(hub, hub$species=='Drosophila melanogaster') @ There are 87 files that match Drosophila melanogaster. If you look at the metadata in hub, you can see that the 7th record representes a GRanges object from UCSC <>= length(hub) head(names(hub)) head(hub$title, n=10) ## then look at a specific slice of the hub object. hub[7] @ So you can retrieve that dm3 file as a \Rcode{GRanges} like this: <>= gr <- hub[[names(hub)[7]]] summary(gr) @ The metadata fields contain the details of file origin and content. <>= metadata(gr) @ Split the \Rclass{GRanges} object by gene name to get a \Rclass{GRangesList} object of transcript ranges grouped by gene. <>= txbygn <- split(gr, gr$name) @ You can now use \Rcode{txbygn} with the \Rcode{summarizeOverlaps} function to prepare a table of read counts for RNA-Seq differential gene expression. Note that before passing \Rcode{txbygn} to \Rcode{summarizeOverlaps}, you should confirm that the seqlevels (chromosome names) in it match those in the BAM file. See \Rcode{?renameSeqlevels}, \Rcode{?keepSeqlevels} and \Rcode{?seqlevels} for examples of renaming seqlevels. %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% \subsection{How to annotate peaks in read coverage} [coming soon...] %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% \subsection{How to prepare a table of read counts for RNA-Seq differential gene expression} Methods for RNA-Seq gene expression analysis generally require a table of counts that summarize the number of reads that overlap or `hit' a particular gene. In this {\it HOWTO} we count with the \Rcode{summarizeOverlaps} function from the \Biocpkg{GenomicAlignments} package and create a count table from the results. Other packages that provide read counting are \Biocpkg{Rsubread} and \Biocpkg{easyRNASeq}. The \Biocpkg{parathyroidSE} package vignette contains a workflow on counting and other common operations required for differential expression analysis. As sample data we use the \Biocpkg{pasillaBamSubset} data package described in the introduction. <>= library(pasillaBamSubset) reads <- c(untrt1=untreated1_chr4(), # single-end reads untrt3=untreated3_chr4()) # paired-end reads @ \Rcode{summarizeOverlaps} requires the name of a BAM file(s) and a {\textit gene model} to count against. See introduction for a quick description of what a \textit{gene models} is. The gene model must match the genome build the reads in the BAM file were aligned to. For the pasilla data this is dm3 Dmelanogaster which is available as a \Bioconductor{} package. Load the package and extract the exon ranges grouped by gene: <>= library(TxDb.Dmelanogaster.UCSC.dm3.ensGene) exbygene <- exonsBy(TxDb.Dmelanogaster.UCSC.dm3.ensGene, "gene") @ \Rcode{exbygene} is a \Rclass{GRangesList} object with one list element per gene in the gene model. \Rcode{summarizeOverlaps} automatically sets a \Rcode{yieldSize} on large BAM files and iterates over them in chunks. When reading paired-end data set the \Rcode{singleEnd} argument to FALSE. See ?\Rfunction{summarizeOverlaps} for details reguarding the count \Rcode{modes} and additional arguments. <>= library(GenomicAlignments) se <- summarizeOverlaps(exbygene, reads, mode="IntersectionNotEmpty") @ The return object is a \Rcode{SummarizedExperiment} with counts accessible with the \Rcode{assays} accessor: <>= class(se) head(table(assays(se)$counts)) @ The count vector is the same length as \Rcode{exbygene}: <>= identical(length(exbygene), length(assays(se)$counts)) @ A copy of \Rcode{exbygene} is stored in the \Rcode{se} object and accessible with \Rcode{rowRanges} accessor: <>= rowRanges(se) @ Two popular packages for RNA-Seq differential gene expression are \Biocpkg{DESeq2} and \Biocpkg{edgeR}. Tables of counts per gene are required for both and can be easily created with a vector of counts. Here we use the counts from our \Rclass{SummarizedExperiment} object: <>= colData(se)$trt <- factor(c("untrt", "untrt"), levels=c("trt", "untrt")) colData(se) library(DESeq2) deseq <- DESeqDataSet(se, design= ~ 1) library(edgeR) edger <- DGEList(assays(se)$counts, group=rownames(colData(se))) @ %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% \subsection{How to summarize junctions from a BAM file containing RNA-Seq reads} As sample data we use the \Biocpkg{pasillaBamSubset} data package described in the introduction. <>= library(pasillaBamSubset) un1 <- untreated1_chr4() # single-end reads library(GenomicAlignments) reads1 <- readGAlignments(un1) reads1 @ For each alignment, the aligner generated a CIGAR string that describes its "geometry", that is, the locations of insertions, deletions and junctions in the alignment. See the SAM Spec available on the SAMtools website for the details (\url{http://samtools.sourceforge.net/}). The \Rcode{summarizeJunctions()} function from the \Biocpkg{GenomicAlignments} package can be used to summarize the junctions in \Rcode{reads1}. <>= junc_summary <- summarizeJunctions(reads1) junc_summary @ See \Rcode{?summarizeJunctions} in the \Biocpkg{GenomicAlignments} package for more information. %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% \subsection{How to get the exon and intron sequences of a given gene} The exon and intron sequences of a gene are essentially the DNA sequences of the introns and exons of all known transcripts of the gene. The first task is to identify all transcripts associated with the gene of interest. Our sample gene is the human TRAK2 which is involved in regulation of endosome-to-lysosome trafficking of membrane cargo. The Entrez gene id is `66008'. <>= trak2 <- "66008" @ The \Biocpkg{TxDb.Hsapiens.UCSC.hg38.knownGene} data package contains the gene model corresponding to the UCSC `Known Genes' track. <>= library(TxDb.Hsapiens.UCSC.hg38.knownGene) txdb <- TxDb.Hsapiens.UCSC.hg38.knownGene @ The transcript ranges for all the genes in the gene model can be extracted with the \Rfunction{transcriptsBy} function from the \Biocpkg{GenomicFeatures} package. They will be returned in a named \Rclass{GRangesList} object containing all the transcripts grouped by gene. In order to keep only the transcripts of the TRAK2 gene we will subset the \Rclass{GRangesList} object using the \Rcode{[[} operator. <>= library(GenomicFeatures) trak2_txs <- transcriptsBy(txdb, by="gene")[[trak2]] trak2_txs @ \Rcode{trak2\_txs} is a \Rclass{GRanges} object with one range per transcript in the TRAK2 gene. The transcript names are stored in the \Rcode{tx\_name} metadata column. We will need them to subset the extracted intron and exon regions: <>= trak2_tx_names <- mcols(trak2_txs)$tx_name trak2_tx_names @ The exon and intron genomic ranges for all the transcripts in the gene model can be extracted with the \Rfunction{exonsBy} and \Rfunction{intronsByTranscript} functions, respectively. Both functions return a \Rclass{GRangesList} object. Then we keep only the exon and intron for the transcripts of the TRAK2 gene by subsetting each \Rclass{GRangesList} object by the TRAK2 transcript names. Extract the exon regions: <>= trak2_exbytx <- exonsBy(txdb, "tx", use.names=TRUE)[trak2_tx_names] elementNROWS(trak2_exbytx) @ ... and the intron regions: <>= trak2_inbytx <- intronsByTranscript(txdb, use.names=TRUE)[trak2_tx_names] elementNROWS(trak2_inbytx) @ Next we want the DNA sequences for these exons and introns. The \Rfunction{getSeq} function from the \Biocpkg{Biostrings} package can be used to query a \Biocpkg{BSgenome} object with a set of genomic ranges and retrieve the corresponding DNA sequences. <>= library(BSgenome.Hsapiens.UCSC.hg38) @ Extract the exon sequences: <>= trak2_ex_seqs <- getSeq(Hsapiens, trak2_exbytx) trak2_ex_seqs trak2_ex_seqs[["uc002uyb.4"]] trak2_ex_seqs[["uc002uyc.2"]] @ ... and the intron sequences: <>= trak2_in_seqs <- getSeq(Hsapiens, trak2_inbytx) trak2_in_seqs trak2_in_seqs[["uc002uyb.4"]] trak2_in_seqs[["uc002uyc.2"]] @ %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% \subsection{How to get the CDS and UTR sequences of genes associated with colorectal cancer} In this {\it HOWTO} we extract the CDS and UTR sequences of genes involved in colorectal cancer. The workflow extends the ideas presented in the previous {\it HOWTO} and suggests an approach for identifying disease-related genes. \subsubsection{Build a gene list} We start with a list of gene or transcript ids. If you do not have pre-defined list one can be created with the \Biocpkg{KEGGREST} and \Biocpkg{KEGGgraph} packages. Create a table of KEGG pathways and ids and search on the term `Colorectal cancer'. <>= library(KEGGREST) li <- keggList("pathway", "hsa") ptag <- names(grep("Colorectal cancer", li, value=TRUE)) ptag tag <- gsub("path:hsa", "", ptag) @ Use the "05210" id to query the KEGG web resource (accesses the currently maintained data). <>= library(KEGGgraph) dest <- tempfile() retrieveKGML(tag, "hsa", dest) @ The suffix of the KEGG id is the Entrez gene id. The \Rfunction{translateKEGGID2GeneID} simply removes the prefix leaving just the Entrez gene ids. <>= crids <- as.character(parseKGML2DataFrame(dest)[,1]) crgenes <- unique(translateKEGGID2GeneID(crids)) head(crgenes) @ \subsubsection{Identify genomic coordinates} The list of gene ids is used to extract genomic positions of the regions of interest. The Known Gene table from UCSC will be the annotation and is available as a \Bioconductor{} package. <>= library(TxDb.Hsapiens.UCSC.hg38.knownGene) txdb <- TxDb.Hsapiens.UCSC.hg38.knownGene @ If an annotation is not available as a \Bioconductor{} annotation package it may be available in \Biocpkg{AnnotationHub}. Additionally, there are functions in \Biocpkg{txdbmaker} which can retrieve data from UCSC and Ensembl to create a \Robject{TxDb}. See \Rcode{?makeTxDbFromUCSC} for more information. As in the previous {\it HOWTO} we need to identify the transcripts corresponding to each gene. The transcript id (or name) is used to isolate the UTR and coding regions of interest. This grouping of transcript by gene is also used to re-group the final sequence results. The \Rcode{transcriptsBy} function outputs both the gene and transcript identifiers which we use to create a map between the two. The \Rcode{map} is a \Robject{CharacterList} with gene ids as names and transcript ids as the list elements. <>= txbygene <- transcriptsBy(txdb, "gene")[crgenes] ## subset on colorectal genes map <- relist(unlist(txbygene, use.names=FALSE)$tx_id, txbygene) map @ Extract the UTR and coding regions. <>= cds <- cdsBy(txdb, "tx") threeUTR <- threeUTRsByTranscript(txdb) fiveUTR <- fiveUTRsByTranscript(txdb) @ Coding and UTR regions may not be present for all transcripts specified in \Rcode{map}. Consequently, the subset results will not be the same length. This length discrepancy must be taken into account when re-listing the final results by gene. <>= txid <- unlist(map, use.names=FALSE) cds <- cds[names(cds) %in% txid] threeUTR <- threeUTR[names(threeUTR) %in% txid] fiveUTR <- fiveUTR[names(fiveUTR) %in% txid] @ Note the different lengths of the subset regions. <>= length(txid) ## all possible transcripts length(cds) length(threeUTR) length(fiveUTR) @ These objects are \Robject{GRangesList}s with the transcript id as the outer list element. <>= cds @ \subsubsection{Extract sequences from BSgenome} The \Rcode{BSgenome} packages contain complete genome sequences for a given organism. Load the \Rcode{BSgenome} package for homo sapiens. <>= library(BSgenome.Hsapiens.UCSC.hg38) genome <- BSgenome.Hsapiens.UCSC.hg38 @ Use \Rfunction{extractTranscriptSeqs} to extract the UTR and coding regions from the \Rcode{BSgenome}. This function retrieves the sequences for an any \Robject{GRanges} or \Robject{GRangesList} (i.e., not just transcripts like the name implies). <>= threeUTR_seqs <- extractTranscriptSeqs(genome, threeUTR) fiveUTR_seqs <- extractTranscriptSeqs(genome, fiveUTR) cds_seqs <- extractTranscriptSeqs(genome, cds) @ The return values are \Robject{DNAStringSet} objects. <>= cds_seqs @ Our final step is to collect the coding and UTR regions (currently organzied by transcript) into groups by gene id. The \Rfunction{relist} function groups the sequences of a \Robject{DNAStringSet} object into a \Robject{DNAStringSetList} object, based on the specified \Rcode{skeleton} argument. The \Rcode{skeleton} must be a list-like object and only its shape (i.e. its element lengths) matters (its exact content is ignored). A simple form of \Rcode{skeleton} is to use a partitioning object that we make by specifying the size of each partition. The partitioning objects are different for each type of region because not all transcripts had a coding or 3' or 5' UTR region defined. <>= lst3 <- relist(threeUTR_seqs, PartitioningByWidth(sum(map %in% names(threeUTR)))) lst5 <- relist(fiveUTR_seqs, PartitioningByWidth(sum(map %in% names(fiveUTR)))) lstc <- relist(cds_seqs, PartitioningByWidth(sum(map %in% names(cds)))) @ There are 239 genes in \Rcode{map} each of which have 1 or more transcripts. The table of element lengths shows how many genes have each number of transcripts. For example, 47 genes have 1 transcript, 48 genes have 2 etc. <>= length(map) table(elementNROWS(map)) @ The lists of DNA sequences all have the same length as \Rcode{map} but one or more of the element lengths may be zero. This would indicate that data were not available for that gene. The tables below show that there was at least 1 coding region available for all genes (i.e., none of the element lengths are 0). However, both the 3' and 5' UTR results have element lengths of 0 which indicates no UTR data were available for that gene. <>= table(elementNROWS(lstc)) table(elementNROWS(lst3)) names(lst3)[elementNROWS(lst3) == 0L] ## genes with no 3' UTR data table(elementNROWS(lst5)) names(lst5)[elementNROWS(lst5) == 0L] ## genes with no 5' UTR data @ %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% \subsection{How to create DNA consensus sequences for read group `families'} The motivation for this {\it HOWTO} comes from a study which explored the dynamics of point mutations. The mutations of interest exist with a range of frequencies in the control group (e.g., 0.1\% - 50\%). PCR and sequencing error rates make it difficult to identify low frequency events (e.g., < 20\%). When a library is prepared with Nextera, random fragments are generated followed by a few rounds of PCR. When the genome is large enough, reads aligning to the same start position are likely descendant from the same template fragment and should have identical sequences. The goal is to elimininate noise by grouping the reads by common start position and discarding those that do not exceed a certain threshold within each family. A new consensus sequence will be created for each read group family. \subsubsection{Sort reads into groups by start position} Load the BAM file into a GAlignments object. <>= library(Rsamtools) bamfile <- system.file("extdata", "ex1.bam", package="Rsamtools") param <- ScanBamParam(what=c("seq", "qual")) library(GenomicAlignments) gal <- readGAlignments(bamfile, use.names=TRUE, param=param) @ Use the \Rfunction{sequenceLayer} function to {\it lay} the query sequences and quality strings on the reference. <>= qseq <- setNames(mcols(gal)$seq, names(gal)) qual <- setNames(mcols(gal)$qual, names(gal)) qseq_on_ref <- sequenceLayer(qseq, cigar(gal), from="query", to="reference") qual_on_ref <- sequenceLayer(qual, cigar(gal), from="query", to="reference") @ Split by chromosome. <>= qseq_on_ref_by_chrom <- splitAsList(qseq_on_ref, seqnames(gal)) qual_on_ref_by_chrom <- splitAsList(qual_on_ref, seqnames(gal)) pos_by_chrom <- splitAsList(start(gal), seqnames(gal)) @ For each chromosome generate one GRanges object that contains unique alignment start positions and attach 3 metadata columns to it: the number of reads, the query sequences, and the quality strings. <>= gr_by_chrom <- lapply(seqlevels(gal), function(seqname) { qseq_on_ref2 <- qseq_on_ref_by_chrom[[seqname]] qual_on_ref2 <- qual_on_ref_by_chrom[[seqname]] pos2 <- pos_by_chrom[[seqname]] qseq_on_ref_per_pos <- split(qseq_on_ref2, pos2) qual_on_ref_per_pos <- split(qual_on_ref2, pos2) nread <- elementNROWS(qseq_on_ref_per_pos) gr_mcols <- DataFrame(nread=unname(nread), qseq_on_ref=unname(qseq_on_ref_per_pos), qual_on_ref=unname(qual_on_ref_per_pos)) gr <- GRanges(Rle(seqname, nrow(gr_mcols)), IRanges(as.integer(names(nread)), width=1)) mcols(gr) <- gr_mcols seqlevels(gr) <- seqlevels(gal) gr }) @ Concatenate all the GRanges objects obtained in (4) together in 1 big GRanges object: <>= gr <- do.call(c, gr_by_chrom) seqinfo(gr) <- seqinfo(gal) @ `gr' is a GRanges object that contains unique alignment start positions: <>= gr[1:6] @ Look at qseq\_on\_ref and qual\_on\_ref. <>= qseq_on_ref qual_on_ref @ 2 reads align to start position 13. Let's have a close look at their sequences: <>= mcols(gr)$qseq_on_ref[[6]] @ and their qualities: <>= mcols(gr)$qual_on_ref[[6]] @ Note that the sequence and quality strings are those projected to the reference so the first letter in those strings are on top of start position 13, the 2nd letter on top of position 14, etc... \subsubsection{Remove low frequency reads} For each start position, remove reads with and under-represented sequence (e.g. threshold = 20\% for the data used here which is low coverage). A unique number is assigned to each unique sequence. This will make future calculations easier and a little bit faster. <>= qseq_on_ref <- mcols(gr)$qseq_on_ref tmp <- unlist(qseq_on_ref, use.names=FALSE) qseq_on_ref_id <- relist(match(tmp, tmp), qseq_on_ref) @ Quick look at `qseq\_on\_ref\_id': It's an IntegerList object with the same length and "shape" as `qseq\_on\_ref'. <>= qseq_on_ref_id @ Remove the under represented ids from each list element of `qseq\_on\_ref\_id': <>= qseq_on_ref_id2 <- endoapply(qseq_on_ref_id, function(ids) ids[countMatches(ids, ids) >= 0.2 * length(ids)]) @ Remove corresponding sequences from `qseq\_on\_ref': <>= tmp <- unlist(qseq_on_ref_id2, use.names=FALSE) qseq_on_ref2 <- relist(unlist(qseq_on_ref, use.names=FALSE)[tmp], qseq_on_ref_id2) @ \subsubsection{Create a consensus sequence for each read group family} Compute 1 consensus matrix per chromosome: <>= split_factor <- rep.int(seqnames(gr), elementNROWS(qseq_on_ref2)) qseq_on_ref2 <- unlist(qseq_on_ref2, use.names=FALSE) qseq_on_ref2_by_chrom <- splitAsList(qseq_on_ref2, split_factor) qseq_pos_by_chrom <- splitAsList(start(gr), split_factor) cm_by_chrom <- lapply(names(qseq_pos_by_chrom), function(seqname) consensusMatrix(qseq_on_ref2_by_chrom[[seqname]], as.prob=TRUE, shift=qseq_pos_by_chrom[[seqname]]-1, width=seqlengths(gr)[[seqname]])) names(cm_by_chrom) <- names(qseq_pos_by_chrom) @ 'cm\_by\_chrom' is a list of consensus matrices. Each matrix has 17 rows (1 per letter in the DNA alphabet) and 1 column per chromosome position. <>= lapply(cm_by_chrom, dim) @ Compute the consensus string from each consensus matrix. We'll put "+" in the strings wherever there is no coverage for that position, and "N" where there is coverage but no consensus. <>= cs_by_chrom <- lapply(cm_by_chrom, function(cm) { ## need to "fix" 'cm' because consensusString() ## doesn't like consensus matrices with columns ## that contain only zeroes (e.g., chromosome ## positions with no coverage) idx <- colSums(cm) == 0L cm["+", idx] <- 1 DNAString(consensusString(cm, ambiguityMap="N")) }) @ The new consensus strings. <>= cs_by_chrom @ %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% \subsection{How to compute binned averages along a genome} In some applications (e.g. visualization), there is the need to compute the average of a variable defined along a genome (a.k.a. genomic variable) for a set of predefined fixed-width regions (sometimes called "bins"). The genomic variable is typically represented as a named \Robject{RleList} object with one list element per chromosome. One such example is coverage. Here we create an artificial genomic variable: <>= library(BSgenome.Scerevisiae.UCSC.sacCer2) set.seed(55) my_var <- RleList( lapply(seqlengths(Scerevisiae), function(seqlen) { tmp <- sample(50L, seqlen, replace=TRUE) %/% 50L Rle(cumsum(tmp - rev(tmp))) } ), compress=FALSE) my_var @ Use the \Rfunction{tileGenome} function to create a set of bins along the genome. <>= bins <- tileGenome(seqinfo(Scerevisiae), tilewidth=100, cut.last.tile.in.chrom=TRUE) @ Compute the binned average for \Rcode{my\_var}: <>= binnedAverage(bins, my_var, "binned_var") @ The bin size can be modified with the \Rcode{tilewidth} argument to \Rfunction{tileGenome}. See \Rcode{?binnedAverage} for additional examples. %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% \section{Session Information} <>= sessionInfo() @ %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% \bibliography{GenomicRanges} \bibliographystyle{plainnat} \end{document}