-
Notifications
You must be signed in to change notification settings - Fork 0
/
GenomicState.Rmd
executable file
·382 lines (297 loc) · 14.2 KB
/
GenomicState.Rmd
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
369
370
371
372
373
374
375
376
377
378
379
380
381
382
---
title: "Introduction to GenomicState"
author:
- name: Leonardo Collado-Torres
affiliation:
- &libd Lieber Institute for Brain Development, Johns Hopkins Medical Campus
- &ccb Center for Computational Biology, Johns Hopkins University
email: lcolladotor@gmail.com
date: "`r doc_date()`"
package: "`r pkg_ver('GenomicState')`"
output:
BiocStyle::html_document
vignette: >
%\VignetteIndexEntry{Introduction to GenomicState}
%\VignetteEngine{knitr::rmarkdown}
%\VignetteEncoding{UTF-8}
---
```{r, include = FALSE}
knitr::opts_chunk$set(
collapse = TRUE,
comment = "#>"
)
```
```{r vignetteSetup, echo=FALSE, message=FALSE, warning = FALSE}
## For links
library("BiocStyle")
## Track time spent on making the vignette
startTime <- Sys.time()
## Bib setup
library("RefManageR")
## Write bibliography information
bib <- c(
R = citation(),
BiocStyle = citation("BiocStyle"),
GenomeInfoDb = RefManageR::BibEntry(
bibtype = "manual",
key = "GenomeInfoDb",
author = "Sonali Arora and Martin Morgan and Marc Carlson and H. Pagès",
title = "GenomeInfoDb: Utilities for manipulating chromosome and other 'seqname' identifiers",
year = 2017, doi = "10.18129/B9.bioc.GenomeInfoDb"
),
knitr = citation("knitr")[3],
GenomicState = citation("GenomicState")[1],
RefManageR = citation("RefManageR")[1],
rmarkdown = citation("rmarkdown")[1],
rtracklayer = citation("rtracklayer")[1],
sessioninfo = citation("sessioninfo"),
testthat = citation("testthat"),
GenomicFeatures = citation("GenomicFeatures"),
bumphunter = citation("bumphunter")[1],
derfinder = citation("derfinder")[1],
AnnotationDbi = citation("AnnotationDbi"),
IRanges = citation("IRanges"),
org.Hs.eg.db = citation("org.Hs.eg.db"),
glue = citation("glue"),
AnnotationHub = citation("AnnotationHub"),
AnnotationHubData = citation("AnnotationHubData"),
GenomicRanges = citation("GenomicRanges")
)
```
# Basics
## Install `r Biocpkg('GenomicState')`
`R` is an open-source statistical environment which can be easily modified to enhance its functionality via packages. `r Biocpkg('GenomicState')` is a `R` package available via Bioconductor. `R` can be installed on any operating system from [CRAN](https://cran.r-project.org/) after which you can install `r Biocpkg('GenomicState')` by using the following commands in your `R` session:
```{r 'install', eval = FALSE}
if (!requireNamespace("BiocManager", quietly = TRUE)) {
install.packages("BiocManager")
}
BiocManager::install("GenomicState")
## Check that you have a valid Bioconductor installation
BiocManager::valid()
```
## Required knowledge
`r Biocpkg('GenomicState')` is based on many other packages and in particular in those that have implemented the infrastructure needed for dealing with annotation data. That is, packages like `r Biocpkg('rtracklayer')` that allow you to import the data. A `r Biocpkg('GenomicState')` user is not expected to deal with those packages directly but will need to be familiar with `r Biocpkg('derfinder')` and `r Biocpkg('derfinderPlot')` to understand the results `r Biocpkg('GenomicState')` generates. Furthermore, it'll be useful for the user to know the syntax of `r Biocpkg('AnnotationHub')` `r Citep(bib[['AnnotationHub']])` in order to query and load the data provided by this package.
If you are asking yourself the question "Where do I start using Bioconductor?" you might be interested in [this blog post](http://lcolladotor.github.io/2014/10/16/startbioc/#.VkOKbq6rRuU).
## Asking for help
As package developers, we try to explain clearly how to use our packages and in which order to use the functions. But `R` and `Bioconductor` have a steep learning curve so it is critical to learn where to ask for help. The blog post quoted above mentions some but we would like to highlight the [Bioconductor support site](https://support.bioconductor.org/) as the main resource for getting help regarding Bioconductor. Other alternatives are available such as creating GitHub issues and tweeting. However, please note that if you want to receive help you should adhere to the [posting guidelines](http://www.bioconductor.org/help/support/posting-guide/). It is particularly critical that you provide a small reproducible example and your session information so package developers can track down the source of the error.
## Citing `r Biocpkg('GenomicState')`
We hope that `r Biocpkg('GenomicState')` will be useful for your research. Please use the following information to cite the package and the overall approach. Thank you!
```{r 'citation'}
## Citation info
citation("GenomicState")
```
# Overview
The `r Biocpkg('GenomicState')` package was developed for speeding up analyses that require these objects and in particular those that rely on Gencode annotation data. The package `r Biocpkg('GenomicState')` provides functions for building `GenomicState` objects from diverse annotation sources such as [`Gencode`](https://www.gencodegenes.org/human/releases.html). It also provides a way to load pre-computed `GenomicState` objects if you are working at [JHPCE](https://jhpce.jhu.edu/). These `GenomicState` objects are normally created using [derfinder::makeGenomicState()](https://rdrr.io/bioc/derfinder/man/makeGenomicState.html) and can be used for annotating regions with [derfinder::annotateRegions()](https://rdrr.io/bioc/derfinder/man/annotateRegions.html) which are in turn used by [derfinderPlot::plotRegionCoverage()](https://rdrr.io/bioc/derfinderPlot/man/plotRegionCoverage.html).
To get started, load the `r Biocpkg('GenomicState')` package.
```{r setup, message = FALSE, warning = FALSE}
library("GenomicState")
```
# AnnotationHub
Using the `GencodeStateHub()` function you can query and download the data from `r Biocpkg('GenomicState')` using `r Biocpkg('AnnotationHub')` `r Citep(bib[['AnnotationHub']])`.
```{r 'annotation_hub'}
## Query AnnotationHub for the GenomicState object for Gencode v31 on
## hg19 coordinates
hub_query_gs_gencode_v31_hg19 <- GenomicStateHub(
version = "31",
genome = "hg19",
filetype = "GenomicState"
)
hub_query_gs_gencode_v31_hg19
## Check the metadata
mcols(hub_query_gs_gencode_v31_hg19)
## Access the file through AnnotationHub
if (length(hub_query_gs_gencode_v31_hg19) == 1) {
hub_gs_gencode_v31_hg19 <- hub_query_gs_gencode_v31_hg19[[1]]
hub_gs_gencode_v31_hg19
}
```
# Using the objects
To show how we can use these objects, first we build those for Gencode version 31 on hg19 coordinates.
```{r 'build_ex_objects'}
## Load the example TxDb object
## or start from scratch with:
## txdb_v31_hg19_chr21 <- gencode_txdb(version = '31', genome = 'hg19',
## chrs = 'chr21')
txdb_v31_hg19_chr21 <- AnnotationDbi::loadDb(
system.file("extdata", "txdb_v31_hg19_chr21.sqlite",
package = "GenomicState"
)
)
## Build the GenomicState and annotated genes
genes_v31_hg19_chr21 <- gencode_annotated_genes(txdb_v31_hg19_chr21)
gs_v31_hg19_chr21 <- gencode_genomic_state(txdb_v31_hg19_chr21)
```
You can alternatively use the files hosted in `r Biocpkg('AnnotationHub')` `r Citep(bib[['AnnotationHub']])` which will be faster in general.
```{r 'obtain_annotation_hub'}
## Create the AnnotationHub object once and re-use it to speed up things
ah <- AnnotationHub::AnnotationHub()
## Find the TxDb object for hg19 Gencode version 31
hub_query_txdb_gencode_v31_hg19 <- GenomicStateHub(
version = "31",
genome = "hg19",
filetype = "TxDb", ah = ah
)
hub_query_txdb_gencode_v31_hg19
## Now the Annotated Genes for hg19 Gencode v31
hub_query_genes_gencode_v31_hg19 <- GenomicStateHub(
version = "31",
genome = "hg19",
filetype = "AnnotatedGenes", ah = ah
)
hub_query_genes_gencode_v31_hg19
## And finally the GenomicState for hg19 Gencode v31
hub_query_gs_gencode_v31_hg19 <- GenomicStateHub(
version = "31",
genome = "hg19",
filetype = "GenomicState", ah = ah
)
hub_query_gs_gencode_v31_hg19
## If you want to access the files use the double bracket AnnotationHub syntax
## to retrieve the R objects from the web.
if (FALSE) {
hub_txdb_gencode_v31_hg19 <- hub_query_txdb_gencode_v31_hg19[[1]]
hub_genes_gencode_v31_hg19 <- hub_query_genes_gencode_v31_hg19[[1]]
hub_gs_gencode_v31_hg19 <- hub_query_gs_gencode_v31_hg19[[1]]
}
```
Next we load a series of related packages that use the objects we created with `r Biocpkg('GenomicState')` or downloaded from `r Biocpkg('AnnotationHub')` `r Citep(bib[['AnnotationHub']])`.
```{r 'load_derfinder', messages = FALSE, warnings = FALSE}
## Load external packages
library("derfinder")
library("derfinderPlot")
library("bumphunter")
library("GenomicRanges")
```
Next we can prepare the needed for running `derfinderPlot::plotRegionCoverage()` where we use the `TxDb` object, the `GenomicState` and the `annotated genes` we prepared for Gencode v31 on hg19.
```{r 'example_plot_prep'}
## Some example regions from derfinder (set the chromosome lengths)
regions <- genomeRegions$regions[1:2]
seqlengths(regions) <- seqlengths(txdb_v31_hg19_chr21)[
names(seqlengths(regions))
]
## Annotate them
nearestAnnotation <- matchGenes(x = regions, subject = genes_v31_hg19_chr21)
annotatedRegions <- annotateRegions(
regions = regions,
genomicState = gs_v31_hg19_chr21$fullGenome, minoverlap = 1
)
## Obtain fullCov object
fullCov <- list("chr21" = genomeDataRaw$coverage)
regionCov <- getRegionCoverage(fullCov = fullCov, regions = regions)
```
And now we can make the example plot as shown below.
```{r 'example_plot'}
## now make the plot
plotRegionCoverage(
regions = regions, regionCoverage = regionCov,
groupInfo = genomeInfo$pop, nearestAnnotation = nearestAnnotation,
annotatedRegions = annotatedRegions, whichRegions = 1:2,
txdb = txdb_v31_hg19_chr21, verbose = FALSE
)
```
# JHPCE
You can also access the data locally using the function `local_metadata()` which works at [JHPCE](https://jhpce.jhu.edu/) or anywhere where you have re-created the files from this package. This returns a `data.frame()` which you can subset. It also inclused the R code for loading the data which you can do using `eval(parse(text = local_metadata()$loadCode))` as shown below.
```{r 'local_metadata'}
## Get the local metadata
meta <- local_metadata()
## Subset to the data of interest, lets say hg19 TxDb for v31
interest <- subset(meta, RDataClass == "TxDb" & Tags == "Gencode:v31:hg19")
## Next you can load the data
if (file.exists(interest$RDataPath)) {
## This only works at JHPCE
eval(parse(text = interest$loadCode))
## Explore the loaded object (would be gencode_v31_hg19_txdb in this case)
gencode_v31_hg19_txdb
}
```
# Build objects
The objects provided by `GenomicState` through `r Biocpkg('AnnotationHub')` `r Citep(bib[['AnnotationHub']])` were built using code like the one included below which is how the Gencode version 23 for hg19 files were built.
```{r 'build_local', eval = FALSE}
outdir <- "gencode"
dir.create(outdir, showWarnings = FALSE)
## Build and save the TxDb object
gencode_v23_hg19_txdb <- gencode_txdb("23", "hg19")
saveDb(gencode_v23_hg19_txdb,
file = file.path(outdir, "gencode_v23_hg19_txdb.sqlite")
)
## Build and save the annotateTranscripts output
gencode_v23_hg19_annotated_genes <- gencode_annotated_genes(
gencode_v23_hg19_txdb
)
save(gencode_v23_hg19_annotated_genes,
file = file.path(outdir, "gencode_v23_hg19_annotated_genes.rda")
)
## Build and save the GenomicState
gencode_v23_hg19_GenomicState <- gencode_genomic_state(
gencode_v23_hg19_txdb
)
save(gencode_v23_hg19_GenomicState,
file = file.path(outdir, "gencode_v23_hg19_GenomicState.rda")
)
```
For more details check the source files:
```{r 'build_source'}
## R commands for building the files:
system.file("scripts", "make-data_gencode_human.R",
package = "GenomicState"
)
## The above file was created by this one:
system.file("scripts", "generate_make_data_gencode_human.R",
package = "GenomicState"
)
```
# Reproducibility
The `r Biocpkg('GenomicState')` package `r Citep(bib[['GenomicState']])` was made possible thanks to:
* R `r Citep(bib[['R']])`
* `r Biocpkg('BiocStyle')` `r Citep(bib[['BiocStyle']])`
* `r Biocpkg('GenomeInfoDb')` `r Citep(bib[['GenomeInfoDb']])`
* `r CRANpkg('knitr')` `r Citep(bib[['knitr']])`
* `r CRANpkg("RefManageR")` `r Citep(bib[["RefManageR"]])`
* `r CRANpkg('rmarkdown')` `r Citep(bib[['rmarkdown']])`
* `r Biocpkg('rtracklayer')` `r Citep(bib[['rtracklayer']])`
* `r CRANpkg('sessioninfo')` `r Citep(bib[['sessioninfo']])`
* `r CRANpkg('testthat')` `r Citep(bib[['testthat']])`
* `r Biocpkg('GenomicFeatures')` `r Citep(bib[['GenomicFeatures']])`
* `r Biocpkg('bumphunter')` `r Citep(bib[['bumphunter']])`
* `r Biocpkg('derfinder')` `r Citep(bib[['derfinder']])`
* `r Biocpkg('AnnotationDbi')` `r Citep(bib[['AnnotationDbi']])`
* `r Biocpkg('IRanges')` `r Citep(bib[['IRanges']])`
* `r Biocpkg('org.Hs.eg.db')` `r Citep(bib[['org.Hs.eg.db']])`
* `r CRANpkg('glue')` `r Citep(bib[['glue']])`
* `r Biocpkg('AnnotationHub')` `r Citep(bib[['AnnotationHub']])`
* `r Biocpkg('AnnotationHubData')` `r Citep(bib[['AnnotationHubData']])`
* `r Biocpkg('GenomicRanges')` `r Citep(bib[['GenomicRanges']])`
Code for creating the vignette
```{r createVignette, eval=FALSE}
## Create the vignette
library("rmarkdown")
system.time(render("GenomicState.Rmd"))
## Extract the R code
library("knitr")
knit("GenomicState.Rmd", tangle = TRUE)
```
Date the vignette was generated.
```{r reproduce1, echo=FALSE}
## Date the vignette was generated
Sys.time()
```
Wallclock time spent generating the vignette.
```{r reproduce2, echo=FALSE}
## Processing time in seconds
totalTime <- diff(c(startTime, Sys.time()))
round(totalTime, digits = 3)
```
`R` session information.
```{r reproduce3, echo=FALSE}
## Session info
library("sessioninfo")
options(width = 120)
session_info()
```
# Bibliography
This vignette was generated using `r Biocpkg('BiocStyle')` `r Citep(bib[['BiocStyle']])`, `r CRANpkg('knitr')` `r Citep(bib[['knitr']])` and `r CRANpkg('rmarkdown')` `r Citep(bib[['rmarkdown']])` running behind the scenes.
Citations made with `r CRANpkg('RefManageR')` `r Citep(bib[['RefManageR']])`.
```{r vignetteBiblio, results = 'asis', echo = FALSE, warning = FALSE, message = FALSE}
## Print bibliography
PrintBibliography(bib, .opts = list(hyperlink = "to.doc", style = "html"))
```