Bioconductor can import diverse sequence-related file types, including fasta, fastq, BAM, VCF, gff, bed, and wig files, among others. Packages support common and advanced sequence manipulation operations such as trimming, transformation, and alignment. Domain-specific analyses include quality assessment, ChIP-seq, differential expression, RNA-seq, and other approaches. Bioconductor includes an interface to the Sequence Read Archive (via the SRAdb package).
This workflow walks through the annotation of a generic set of ranges with Bioconductor packages. The ranges can be any user-defined region of interest or can be from a public file.
As a first step, data are put into a GRanges object so we can take advantage of overlap operations and store identifiers as metadata columns.
The first set of ranges are variants from a dbSNP Variant Call Format (VCF) file. This file can be downloaded from the ftp site at NCBI ftp://ftp.ncbi.nlm.nih.gov/snp/ and imported with readVcf() from the VariantAnnotation package. Alternatively, the file is available as a pre-parsed VCF object in the AnnotationHub.
library(VariantAnnotation)
library(AnnotationHub)
The Hub returns a VcfFile object with a reference to the file on disk.
hub <- AnnotationHub()
## updating AnnotationHub metadata: retrieving 1 resource
## snapshotDate(): 2015-08-26
fl <- hub[['AH47004']]
## loading from cache '/var/lib/jenkins/.AnnotationHub/52446'
## '/var/lib/jenkins/.AnnotationHub/52447'
fl
## class: VcfFile
## path: /var/lib/jenkins/.AnnotationHub/52446
## index: /var/lib/jenkins/.AnnotationHub/52447
## isOpen: FALSE
## yieldSize: NA
Read the data into a VCF object:
vcf <- readVcf(fl, "hg19")
dim(vcf)
## [1] 100242 0
Overlap operations require that seqlevels and the genome of the objects match. Here we modify the VCF seqlevels to match the TxDb.
library(TxDb.Hsapiens.UCSC.hg19.knownGene)
txdb_hg19 <- TxDb.Hsapiens.UCSC.hg19.knownGene
head(seqlevels(txdb_hg19))
## [1] "chr1" "chr2" "chr3" "chr4" "chr5" "chr6"
seqlevels(vcf)
## [1] "1" "2" "3" "4" "5" "6" "7" "8" "9" "10" "11" "12" "13" "14"
## [15] "15" "16" "17" "18" "19" "20" "21" "22" "X" "Y" "MT"
seqlevels(vcf) <- paste0("chr", seqlevels(vcf))
Sanity check to confirm we have matching seqlevels.
intersect(seqlevels(txdb_hg19), seqlevels(vcf))
## [1] "chr1" "chr2" "chr3" "chr4" "chr5" "chr6" "chr7" "chr8"
## [9] "chr9" "chr10" "chr11" "chr12" "chr13" "chr14" "chr15" "chr16"
## [17] "chr17" "chr18" "chr19" "chr20" "chr21" "chr22" "chrX" "chrY"
The genomes already match so no change is needed.
unique(genome(txdb_hg19))
## [1] "hg19"
unique(genome(vcf))
## [1] "hg19"
We are only interested in the standard chromosomes so we drop the rest.
txdb_hg19 <- keepStandardChromosomes(txdb_hg19)
vcf <- keepStandardChromosomes(vcf)
The GRanges in a VCF object is extracted with ‘rowRanges()’.
gr_hg19 <- rowRanges(vcf)
The second set of ranges is a user-defined region of chromosome 4 in mouse. The idea here is that any region, known or unknown, can be annotated with the following steps.
Load the TxDb package and keep only the standard chromosomes.
library(TxDb.Mmusculus.UCSC.mm10.ensGene)
txdb_mm10 <- keepStandardChromosomes(TxDb.Mmusculus.UCSC.mm10.ensGene)
We are creating the GRanges from scratch and can specify the seqlevels (chromosome names) to match the TxDb.
head(seqlevels(txdb_mm10))
## [1] "chr1" "chr2" "chr3" "chr4" "chr5" "chr6"
gr_mm10 <- GRanges("chr4", IRanges(c(4000000, 107889000), width=1000))
Now assign the genome.
unique(genome(txdb_mm10))
## [1] "mm10"
genome(gr_mm10) <- "mm10"
locateVariants() in the VariantAnnotation package annotates ranges with transcript, exon, cds and gene ID’s from a TxDb. Various extractions are performed on the TxDb (exonsBy(), transcripts(), cdsBy(), etc.) and the result is overlapped with the ranges. An appropriate GRangesList can also be supplied as the annotation. Different variants such as ‘coding’, ‘fiveUTR’, ‘threeUTR’, ‘spliceSite’, ‘intron’, ‘promoter’, and ‘intergenic’ can be searched for by passing the appropriate constructor as the ‘region’ argument. See ?locateVariants for details.
loc_hg19 <- locateVariants(gr_hg19, txdb_hg19, AllVariants())
## 'select()' returned many:1 mapping between keys and columns
## 'select()' returned many:1 mapping between keys and columns
## 'select()' returned many:1 mapping between keys and columns
## 'select()' returned many:1 mapping between keys and columns
## 'select()' returned many:1 mapping between keys and columns
## 'select()' returned many:1 mapping between keys and columns
table(loc_hg19$LOCATION)
##
## spliceSite intron fiveUTR threeUTR coding intergenic
## 326 211920 1249 5083 11523 36207
## promoter
## 12343
loc_mm10 <- locateVariants(gr_mm10, txdb_mm10, AllVariants())
## 'select()' returned 1:1 mapping between keys and columns
## 'select()' returned 1:1 mapping between keys and columns
## 'select()' returned 1:1 mapping between keys and columns
table(loc_mm10$LOCATION)
##
## spliceSite intron fiveUTR threeUTR coding intergenic
## 6 1 0 0 0 0
## promoter
## 12
The ID’s returned from locateVariants() can be used in select() to map to ID’s in other annotation packages.
library(org.Hs.eg.db)
cols <- c("UNIPROT", "PFAM")
keys <- na.omit(unique(loc_hg19$GENEID))
head(select(org.Hs.eg.db, keys, cols, keytype="ENTREZID"))
## ENTREZID UNIPROT PFAM
## 1 54991 A0A024R082 PF14946
## 2 54991 Q96HA4 PF14946
## 3 54991 Q5T2W9 PF14946
## 4 116983 Q8WTZ1 PF00169
## 5 116983 Q8WTZ1 PF01412
## 6 116983 Q8WTZ1 PF12796
The ‘keytype’ argument specifies that the mouse TxDb contains Ensembl instead of Entrez gene id’s.
library(org.Mm.eg.db)
keys <- unique(loc_mm10$GENEID)
head(select(org.Mm.eg.db, keys, cols, keytype="ENSEMBL"))
## ENSEMBL UNIPROT PFAM
## 1 ENSMUSG00000028236 Q7TQA3 PF00106
## 2 ENSMUSG00000028608 Q8BHG2 PF05907
Files stored in the AnnotationHub have been pre-processed into ranged-based R objects such as a GRanges, GAlignments and VCF. The positions in our GRanges can be overlapped with the ranges in the AnnotationHub files. This allows for easy subsetting of multiple files, resulting in only the ranges of interest.
Create a ‘hub’ from AnnotationHub and filter the files based on organism and genome build.
hub <- AnnotationHub()
## snapshotDate(): 2015-08-26
hub_hg19 <- subset(hub,
(hub$species == "Homo sapiens") & (hub$genome == "hg19"))
length(hub_hg19)
## [1] 23787
Iterate over the first 3 files and extract ranges that overlap ‘gr_hg19’.
ov_hg19 <- lapply(1:3, function(i) subsetByOverlaps(hub_hg19[[i]], gr_hg19))
Inspect the results.
names(ov_hg19) <- names(hub_hg19)[1:3]
lapply(ov_hg19, head, n=3)
## $AH3166
## GRanges object with 3 ranges and 5 metadata columns:
## seqnames ranges strand |
## <Rle> <IRanges> <Rle> |
## [1] chr14 [ 23388231, 23388425] - |
## [2] chr11 [118436595, 118436793] - |
## [3] chr9 [ 21994302, 21994477] - |
## name score level
## <character> <integer> <numeric>
## [1] chr14:23388231:23388425:-:1.08:0.993703 0 370.1525
## [2] chr11:118436595:118436793:-:2.21:0.999420 0 255.1222
## [3] chr9:21994302:21994477:-:0.14:1.000000 0 156.7853
## signif score2
## <numeric> <integer>
## [1] 5.15e-09 0
## [2] 8.41e-09 0
## [3] 1.17e-08 0
## -------
## seqinfo: 24 sequences (1 circular) from hg19 genome
##
## $AH3912
## GRanges object with 3 ranges and 5 metadata columns:
## seqnames ranges strand | name score
## <Rle> <IRanges> <Rle> | <character> <integer>
## [1] chr1 [1048597, 1049288] * | . 0
## [2] chr1 [1050193, 1050947] * | . 0
## [3] chr1 [1051101, 1053274] * | . 0
## signalValue pValue qValue
## <numeric> <numeric> <numeric>
## [1] 12.5496 21.9043 -1
## [2] 10.1759 15.7118 -1
## [3] 897.2450 324.0000 -1
## -------
## seqinfo: 23 sequences from hg19 genome
##
## $AH3913
## GRanges object with 3 ranges and 6 metadata columns:
## seqnames ranges strand | name score
## <Rle> <IRanges> <Rle> | <character> <integer>
## [1] chr1 [ 1048860, 1049010] * | . 0
## [2] chr1 [ 9263720, 9263870] * | . 0
## [3] chr1 [17027700, 17027850] * | . 0
## signalValue pValue qValue peak
## <numeric> <numeric> <numeric> <integer>
## [1] 60 21.9043 -1 -1
## [2] 238 210.8220 -1 -1
## [3] 102 48.3008 -1 -1
## -------
## seqinfo: 23 sequences from hg19 genome
Annotating the mouse ranges in the same fashion is left as an exercise.
For the set of dbSNP variants that fall in coding regions, amino acid changes can be computed. The output contains one line for each variant-transcript match which can result in multiple lines for each variant.
library(BSgenome.Hsapiens.UCSC.hg19)
head(predictCoding(vcf, txdb_hg19, Hsapiens), 3)
## GRanges object with 3 ranges and 17 metadata columns:
## seqnames ranges strand | paramRangeID
## <Rle> <IRanges> <Rle> | <factor>
## rs397514721 chr1 [1233203, 1233203] - | <NA>
## rs397514721 chr1 [1233203, 1233203] - | <NA>
## rs397514721 chr1 [1233203, 1233203] - | <NA>
## REF ALT QUAL FILTER
## <DNAStringSet> <DNAStringSetList> <numeric> <character>
## rs397514721 T A <NA> .
## rs397514721 T A <NA> .
## rs397514721 T A <NA> .
## varAllele CDSLOC PROTEINLOC QUERYID
## <DNAStringSet> <IRanges> <IntegerList> <integer>
## rs397514721 T [ 317, 317] 106 55
## rs397514721 T [1001, 1001] 334 55
## rs397514721 T [1127, 1127] 376 55
## TXID CDSID GENEID CONSEQUENCE
## <character> <IntegerList> <character> <factor>
## rs397514721 4150 12185 116983 nonsynonymous
## rs397514721 4151 12185 116983 nonsynonymous
## rs397514721 4152 12185 116983 nonsynonymous
## REFCODON VARCODON REFAA VARAA
## <DNAStringSet> <DNAStringSet> <AAStringSet> <AAStringSet>
## rs397514721 GAG GTG E V
## rs397514721 GAG GTG E V
## rs397514721 GAG GTG E V
## -------
## seqinfo: 24 sequences from hg19 genome; no seqlengths
sessionInfo()
## R version 3.2.2 (2015-08-14)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu precise (12.04.4 LTS)
##
## locale:
## [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
## [3] LC_TIME=en_US.UTF-8 LC_COLLATE=C
## [5] LC_MONETARY=en_US.UTF-8 LC_MESSAGES=en_US.UTF-8
## [7] LC_PAPER=en_US.UTF-8 LC_NAME=C
## [9] LC_ADDRESS=C LC_TELEPHONE=C
## [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C
##
## attached base packages:
## [1] stats4 parallel stats graphics grDevices utils datasets
## [8] methods base
##
## other attached packages:
## [1] BSgenome.Hsapiens.UCSC.hg19_1.4.0
## [2] BSgenome_1.38.0
## [3] rtracklayer_1.30.1
## [4] org.Mm.eg.db_3.2.3
## [5] org.Hs.eg.db_3.2.3
## [6] RSQLite_1.0.0
## [7] DBI_0.3.1
## [8] TxDb.Mmusculus.UCSC.mm10.ensGene_3.2.2
## [9] TxDb.Hsapiens.UCSC.hg19.knownGene_3.2.2
## [10] GenomicFeatures_1.22.4
## [11] AnnotationDbi_1.32.0
## [12] AnnotationHub_2.2.1
## [13] VariantAnnotation_1.16.3
## [14] Rsamtools_1.22.0
## [15] Biostrings_2.38.0
## [16] XVector_0.10.0
## [17] SummarizedExperiment_1.0.0
## [18] Biobase_2.30.0
## [19] GenomicRanges_1.22.0
## [20] GenomeInfoDb_1.6.1
## [21] IRanges_2.4.1
## [22] S4Vectors_0.8.1
## [23] BiocGenerics_0.16.0
##
## loaded via a namespace (and not attached):
## [1] Rcpp_0.12.1 BiocInstaller_1.20.0
## [3] formatR_1.2.1 futile.logger_1.4.1
## [5] bitops_1.0-6 futile.options_1.0.0
## [7] tools_3.2.2 zlibbioc_1.16.0
## [9] biomaRt_2.26.0 digest_0.6.8
## [11] evaluate_0.8 shiny_0.12.2
## [13] curl_0.9.3 yaml_2.1.13
## [15] stringr_1.0.0 httr_1.0.0
## [17] knitr_1.11 R6_2.1.1
## [19] XML_3.98-1.3 BiocParallel_1.4.0
## [21] rmarkdown_0.8.1 lambda.r_1.1.7
## [23] magrittr_1.5 htmltools_0.2.6
## [25] GenomicAlignments_1.6.1 xtable_1.8-0
## [27] mime_0.4 interactiveDisplayBase_1.8.0
## [29] httpuv_1.3.3 stringi_1.0-1
## [31] RCurl_1.95-4.7
Exercise 1: VCF header and reading data subsets.
VCF files can be large and it’s often the case that only a subset of variables or genomic positions are of interest. The scanVcfHeader() function in the VariantAnnotation package retrieves header information from a VCF file. Based on the information returned from scanVcfHeader() a ScanVcfParam() object can be created to read in a subset of data from a VCF file. * Use scanVcfHeader() to inspect the header information in the ‘chr22.vcf.gz’ file in VariantAnnotation package. * Select a few ‘info’ or ‘geno’ variables and create a ScanVcfParam object. * Use the ScanVcfParam object as the ‘param’ argument to readVcf() to read in a subset of data. Note that the header() accessor operates on VCF objects in the R workspace. Try header(vcf) on the dbSNP file from AnnotationHub.
Exercise 2: Annotate the mouse ranges in ‘gr_mm10’ with AnnotationHub files. * Create a new ‘hub’ and filter on organism. * Isolate the files for the appropriate genome build and perform overlaps.
Exercise 3: Annotate a gene range from Saccharomyces Scerevisiae. * Load TxDb.Scerevisiae.UCSC.sacCer3.sgdGene and extract the gene ranges. (Hint: use transcriptsBy() and range()). * Isolate the range for gene “YBL086C”. * Create a new ‘hub’ from AnnotationHub and filter by organism. (You should see >= 39 files.) * Select the files for ‘sacCer3’ and perform overlaps.
[ Back to top ]