1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368
|
---
title: "AnnotationHub How-To's"
output:
BiocStyle::html_document:
toc: true
vignette: >
% \VignetteIndexEntry{AnnotationHub: AnnotationHub HOW TO's}
% \VignetteDepends{AnnotationHub, txdbmaker, Rsamtools}
% \VignetteEngine{knitr::rmarkdown}
---
```{r style, echo = FALSE, results = 'asis', warning=FALSE}
options(width=100)
suppressPackageStartupMessages({
## load here to avoid noise in the body of the vignette
library(AnnotationHub)
library(txdbmaker)
library(Rsamtools)
library(VariantAnnotation)
})
BiocStyle::markdown()
```
**Package**: `r Biocpkg("AnnotationHub")`<br />
**Authors**: `r packageDescription("AnnotationHub")[["Author"]] `<br />
**Modified**: Mon March 18 2024<br />
**Compiled**: `r date()`
# Accessing Genome-Scale Data
## Non-model organism gene annotations
_Bioconductor_ offers pre-built `org.*` annotation packages for model
organisms, with their use described in the
[OrgDb](http://bioconductor.org/help/workflows/annotation/Annotation_Resources/#OrgDb)
section of the Annotation work flow. Here we discover available `OrgDb`
objects for less-model organisms
```{r less-model-org}
library(AnnotationHub)
ah <- AnnotationHub()
query(ah, "OrgDb")
orgdb <- query(ah, c("OrgDb", "maintainer@bioconductor.org"))[[1]]
```
The object returned by AnnotationHub is directly usable with the
`select()` interface, e.g., to discover the available keytypes for
querying the object, the columns that these keytypes can map to, and
finally selecting the SYMBOL and GENENAME corresponding to the first 6
ENTREZIDs
```{r less-model-org-select}
library(AnnotationDbi)
AnnotationDbi::keytypes(orgdb)
AnnotationDbi::columns(orgdb)
egid <- head(keys(orgdb, "ENTREZID"))
AnnotationDbi::select(orgdb, egid, c("SYMBOL", "GENENAME"), "ENTREZID")
```
## Roadmap Epigenomics Project
All Roadmap Epigenomics files are hosted
[here](http://egg2.wustl.edu/roadmap/data/byFileType/). If one had to
download these files on their own, one would navigate through the web
interface to find useful files, then use something like the following
_R_ code.
```{r, eval=FALSE}
url <- "http://egg2.wustl.edu/roadmap/data/byFileType/peaks/consolidated/broadPeak/E001-H3K4me1.broadPeak.gz"
filename <- basename(url)
download.file(url, destfile=filename)
if (file.exists(filename))
data <- import(filename, format="bed")
```
This would have to be repeated for all files, and the onus would lie
on the user to identify, download, import, and manage the local disk
location of these files.
`r Biocpkg("AnnotationHub")` reduces this task to just a few lines of _R_ code
```{r results='hide'}
library(AnnotationHub)
ah = AnnotationHub()
epiFiles <- query(ah, "EpigenomeRoadMap")
```
A look at the value returned by `epiFiles` shows us that
`r length(epiFiles)` roadmap resources are available via
`r Biocpkg("AnnotationHub")`. Additional information about
the files is also available, e.g., where the files came from
(dataprovider), genome, species, sourceurl, sourcetypes.
```{r}
epiFiles
```
A good sanity check to ensure that we have files only from the Roadmap Epigenomics
project is to check that all the files in the returned smaller hub object
come from _Homo sapiens_ and the `r unique(epiFiles$genome)` genome
```{r}
unique(epiFiles$species)
unique(epiFiles$genome)
```
Broadly, one can get an idea of the different files from this project
looking at the sourcetype
```{r}
table(epiFiles$sourcetype)
```
To get a more descriptive idea of these different files one can use:
```{r}
sort(table(epiFiles$description), decreasing=TRUE)
```
The 'metadata' provided by the Roadmap Epigenomics Project is also
available. Note that the information displayed about a hub with a
single resource is quite different from the information displayed when
the hub references more than one resource.
```{r}
metadata.tab <- query(ah , c("EpigenomeRoadMap", "Metadata"))
metadata.tab
```
So far we have been exploring information about resources, without
downloading the resource to a local cache and importing it into R.
One can retrieve the resource using `[[` as indicated at the
end of the show method
```{r echo=FALSE, results='hide'}
metadata.tab <- ah[["AH41830"]]
```
```{r}
metadata.tab <- ah[["AH41830"]]
```
The metadata.tab file is returned as a _data.frame_. The first 6 rows
of the first 5 columns are shown here:
```{r}
metadata.tab[1:6, 1:5]
```
One can keep constructing different queries using multiple arguments to
trim down these `r length(epiFiles)` to get the files one wants.
For example, to get the ChIP-Seq files for consolidated epigenomes,
one could use
```{r}
bpChipEpi <- query(ah , c("EpigenomeRoadMap", "broadPeak", "chip", "consolidated"))
```
To get all the bigWig signal files, one can query the hub using
```{r}
allBigWigFiles <- query(ah, c("EpigenomeRoadMap", "BigWig"))
```
To access the 15 state chromatin segmentations, one can use
```{r}
seg <- query(ah, c("EpigenomeRoadMap", "segmentations"))
```
If one is interested in getting all the files related to one sample
```{r}
E126 <- query(ah , c("EpigenomeRoadMap", "E126", "H3K4ME2"))
E126
```
Hub resources can also be selected using `$`, `subset()`, and
`BiocHubsShiny()`; see the main
[_AnnotationHub_ vignette](AnnotationHub.html) for additional detail.
Hub resources are imported as the appropriate _Bioconductor_ object
for use in further analysis. For example, peak files are returned as
_GRanges_ objects.
```{r echo=FALSE, results='hide'}
peaks <- E126[['AH29817']]
```
```{r}
peaks <- E126[['AH29817']]
seqinfo(peaks)
```
BigWig files are returned as _BigWigFile_ objects. A _BigWigFile_ is a
reference to a file on disk; the data in the file can be read in using
`rtracklayer::import()`, perhaps querying these large files for
particular genomic regions of interest as described on the help page
`?import.bw`.
Each record inside `r Biocpkg("AnnotationHub")` is associated with a
unique identifier. Most _GRanges_ objects returned by
`r Biocpkg("AnnotationHub")` contain the unique AnnotationHub identifier of
the resource from which the _GRanges_ is derived. This can come handy
when working with the _GRanges_ object for a while, and additional
information about the object (e.g., the name of the file in the cache,
or the original sourceurl for the data underlying the resource) that
is being worked with.
```{r}
metadata(peaks)
ah[metadata(peaks)$AnnotationHubName]$sourceurl
```
## Ensembl GTF and FASTA files for TxDb gene models and sequence queries
_Bioconductor_ represents gene models using 'transcript'
databases. These are available via packages such as
`r Biocannopkg("TxDb.Hsapiens.UCSC.hg38.knownGene")`
or can be constructed using functions such as
`r Biocpkg("txdbmaker")`::`makeTxDbFromBiomart()`.
_AnnotationHub_ provides an easy way to work with gene models
published by Ensembl. Let's see what Ensembl's Release-94 has in terms
of data for pufferfish, _Takifugu rubripes_.
```{r takifugu-gene-models}
query(ah, c("Takifugu", "release-94"))
```
We see that there is a GTF file descrbing gene models, as well as
various DNA sequences. Let's retrieve the GTF and top-level DNA
sequence files. The GTF file is imported as a _GRanges_ instance, the
DNA sequence as a twobit file.
```{r takifugu-data}
gtf <- ah[["AH64858"]]
dna <- ah[["AH66116"]]
head(gtf, 3)
dna
head(seqlevels(dna))
```
Let's identify the 25 longest DNA sequences, and keep just the
annotations on these scaffolds.
```{r takifugu-seqlengths}
keep <- names(tail(sort(seqlengths(dna)), 25))
gtf_subset <- gtf[seqnames(gtf) %in% keep]
```
It is trivial to make a TxDb instance of this subset (or of the entire
gtf)
```{r takifugu-txdb}
library(txdbmaker) # for makeTxDbFromGRanges
txdb <- makeTxDbFromGRanges(gtf_subset)
```
and to use that in conjunction with the DNA sequences, e.g., to find
exon sequences of all annotated genes.
```{r takifugu-exons}
library(Rsamtools) # for getSeq,FaFile-method
exons <- exons(txdb)
length(exons)
getSeq(dna, exons)
```
There is a one-to-one mapping between the genomic ranges contained in
`exons` and the DNA sequences returned by `getSeq()`.
Some difficulties arise when working with this partly assembled genome
that require more advanced GenomicRanges skills, see the
`r Biocpkg("GenomicRanges")` vignettes, especially "_GenomicRanges_
HOWTOs" and "An Introduction to _GenomicRanges_".
## liftOver to map between genome builds
Suppose we wanted to lift features from one genome build to another,
e.g., because annotations were generated for hg19 but our experimental
analysis used hg18. We know that UCSC provides 'liftover' files for
mapping between genome builds.
In this example, we will take our broad Peak _GRanges_ from E126 which
comes from the 'hg19' genome, and lift over these features to their
'hg38' coordinates.
```{r}
chainfiles <- query(ah , c("hg38", "hg19", "chainfile"))
chainfiles
```
We are interested in the file that lifts over features from hg19 to
hg38 so lets download that using
```{r echo=FALSE, results='hide'}
chain <- chainfiles[['AH14150']]
```
```{r}
chain <- chainfiles[['AH14150']]
chain
```
Perform the liftOver operation using `rtracklayer::liftOver()`:
```{r}
library(rtracklayer)
gr38 <- liftOver(peaks, chain)
```
This returns a _GRangeslist_; update the genome of the result to get
the final result
```{r}
genome(gr38) <- "hg38"
gr38
```
## Working with dbSNP Variants
One may also be interested in working with common germline variants with
evidence of medical interest. This information is available at
[NCBI](https://www.ncbi.nlm.nih.gov/variation/docs/human_variation_vcf/).
Query the dbDNP files in the hub:
```{r echo=FALSE, results='hide', message=FALSE}
query(ah, c("GRCh38", "dbSNP", "VCF" ))
vcf <- ah[['AH57960']]
```
This returns a _VcfFile_ which can be read in using `r
Biocpkg("VariantAnnotation")`; because VCF files can be large, `readVcf()`
supports several strategies for importing only relevant parts of the file
(e.g., particular genomic locations, particular features of the variants), see
`?readVcf` for additional information.
```{r message=FALSE}
variants <- readVcf(vcf, genome="hg19")
variants
```
`rowRanges()` returns information from the CHROM, POS and ID fields of the VCF
file, represented as a _GRanges_ instance
```{r}
rowRanges(variants)
```
Note that the broadPeaks files follow the UCSC chromosome naming convention,
and the vcf data follows the NCBI style of chromosome naming convention.
To bring these ranges in the same chromosome
naming convention (ie UCSC), we would use
```{r}
seqlevelsStyle(variants) <-seqlevelsStyle(peaks)
```
And then finally to find which variants overlap these broadPeaks we would use:
```{r}
overlap <- findOverlaps(variants, peaks)
overlap
```
Some insight into how these results can be interpretted comes from
looking a particular peak, e.g., the 3852nd peak
```{r}
idx <- subjectHits(overlap) == 3852
overlap[idx]
```
There are three variants overlapping this peak; the coordinates of the
peak and the overlapping variants are
```{r}
peaks[3852]
rowRanges(variants)[queryHits(overlap[idx])]
```
# sessionInfo
```{r}
sessionInfo()
```
|