This page describes how to construct a GeneTrack/FeatureTrack from TxDb or GRanges.

suppressPackageStartupMessages({
    library(TnT)
    library(GenomicFeatures)
    library(TxDb.Hsapiens.UCSC.hg19.knownGene)
    library(EnsDb.Hsapiens.v75)
})

Gene Track From TxDb

A GeneTrack can be easily constructed from a TxDb.

txdb <- TxDb.Hsapiens.UCSC.hg19.knownGene
gtrack <- TnT::GeneTrackFromTxDb(txdb)
gtrack
## A GeneTrack 
## | Label: txdb
## | Background:    missing, use 'white'
## | Height:    100
## | Data:  GeneTrackData object with 23056 ranges and 5 metadata columns:
## |          seqnames                 ranges strand |      tooltip       color
## |             <Rle>              <IRanges>  <Rle> | <data.frame> <character>
## |        1    chr19 [ 58858172,  58874214]      - |            1       black
## |       10     chr8 [ 18248755,  18258723]      + |           10       black
## |      100    chr20 [ 43248163,  43280376]      - |          100       black
## |     1000    chr18 [ 25530930,  25757445]      - |         1000       black
## |    10000     chr1 [243651535, 244006886]      - |        10000       black
## |      ...      ...                    ...    ... .          ...         ...
## |     9991     chr9 [114979995, 115095944]      - |         9991       black
## |     9992    chr21 [ 35736323,  35743440]      + |         9992       black
## |     9993    chr22 [ 19023795,  19109967]      - |         9993       black
## |     9994     chr6 [ 90539619,  90584155]      + |         9994       black
## |     9997    chr22 [ 50961997,  50964905]      - |         9997       black
## |                key display_label          id
## |          <integer>   <character> <character>
## |        1         1           < 1           1
## |       10         2          10 >          10
## |      100         3         < 100         100
## |     1000         4        < 1000        1000
## |    10000         5       < 10000       10000
## |      ...       ...           ...         ...
## |     9991     23052        < 9991        9991
## |     9992     23053        9992 >        9992
## |     9993     23054        < 9993        9993
## |     9994     23055        9994 >        9994
## |     9997     23056        < 9997        9997
## |    -------
## |    seqinfo: 93 sequences (1 circular) from hg19 genome

To show the track, put the track into a TnTBoard or TnTGenome, optionally specifing the view range with a GRanges object.

TnTGenome(gtrack, view.range = GRanges("chr13", IRanges(32889617, 32973809)) * 0.7)
## Warning in readLines(con): incomplete final line found on '/home/jialin/R/
## x86_64-suse-linux-gnu-library/3.4/TnT/htmlwidgets/trackWidget.yaml'



Feature Track

Feature Track Showing Genes

FeatureTrack is an alias of GeneTrack, in the sense that this display method can be applied to general overlapping genomic features.

Function FeatureTrack can construct a GeneTrack/FeatureTrack from a GRanges object, and gives you more control of the display (i.e. color, feature labels, tooltip).

For example, we first extract the genes as a GRanges object from a EnsDb.

gene <- genes(EnsDb.Hsapiens.v75)
head(gene)
## GRanges object with 6 ranges and 6 metadata columns:
##                   seqnames         ranges strand |         gene_id
##                      <Rle>      <IRanges>  <Rle> |     <character>
##   ENSG00000223972        1 [11869, 14412]      + | ENSG00000223972
##   ENSG00000227232        1 [14363, 29806]      - | ENSG00000227232
##   ENSG00000243485        1 [29554, 31109]      + | ENSG00000243485
##   ENSG00000237613        1 [34554, 36081]      - | ENSG00000237613
##   ENSG00000268020        1 [52473, 54936]      + | ENSG00000268020
##   ENSG00000240361        1 [62948, 63887]      + | ENSG00000240361
##                     gene_name gene_biotype seq_coord_system      symbol
##                   <character>  <character>      <character> <character>
##   ENSG00000223972     DDX11L1   pseudogene       chromosome     DDX11L1
##   ENSG00000227232      WASH7P   pseudogene       chromosome      WASH7P
##   ENSG00000243485  MIR1302-10      lincRNA       chromosome  MIR1302-10
##   ENSG00000237613     FAM138A      lincRNA       chromosome     FAM138A
##   ENSG00000268020      OR4G4P   pseudogene       chromosome      OR4G4P
##   ENSG00000240361     OR4G11P   pseudogene       chromosome     OR4G11P
##                                            entrezid
##                                              <list>
##   ENSG00000223972               100287596,100287102
##   ENSG00000227232                  100287171,653635
##   ENSG00000243485 100422919,100422834,100422831,...
##   ENSG00000237613              654835,645520,641702
##   ENSG00000268020                                NA
##   ENSG00000240361                                NA
##   -------
##   seqinfo: 273 sequences from GRCh37 genome

Then we construct the GeneTrack using this GRanges object with FeatureTrack, specifing our custom feature labels, feature colors and tooltips. Finally show the track with TnTGenome.

ensGeneTrack <- TnT::FeatureTrack(gene, tooltip = as.data.frame(gene),
                       names = paste(gene$symbol, " (", gene$gene_biotype, ")", sep = ""),
                       color = TnT::mapcol(gene$gene_biotype, palette.fun = grDevices::rainbow))
TnTGenome(ensGeneTrack, view.range = gene[gene$symbol == "BRCA2"][1] * .7)
## Warning in readLines(con): incomplete final line found on '/home/jialin/R/
## x86_64-suse-linux-gnu-library/3.4/TnT/htmlwidgets/trackWidget.yaml'



Feature Track Showing Alignment

TODO



Session Info

sessionInfo()
## R version 3.4.1 (2017-06-30)
## Platform: x86_64-suse-linux-gnu (64-bit)
## Running under: openSUSE Tumbleweed
## 
## Matrix products: default
## BLAS: /usr/lib64/R/lib/libRblas.so
## LAPACK: /usr/lib64/R/lib/libRlapack.so
## 
## locale:
##  [1] LC_CTYPE=en_US.UTF-8       LC_NUMERIC=C              
##  [3] LC_TIME=en_US.UTF-8        LC_COLLATE=en_US.UTF-8    
##  [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] parallel  stats4    methods   stats     graphics  grDevices utils    
## [8] datasets  base     
## 
## other attached packages:
##  [1] EnsDb.Hsapiens.v75_2.99.0              
##  [2] ensembldb_2.1.13                       
##  [3] AnnotationFilter_1.1.9                 
##  [4] TxDb.Hsapiens.UCSC.hg19.knownGene_3.2.2
##  [5] GenomicFeatures_1.30.3                 
##  [6] AnnotationDbi_1.40.0                   
##  [7] Biobase_2.38.0                         
##  [8] TnT_1.3.1                              
##  [9] GenomicRanges_1.30.3                   
## [10] GenomeInfoDb_1.14.0                    
## [11] IRanges_2.12.0                         
## [12] S4Vectors_0.16.0                       
## [13] BiocGenerics_0.24.0                    
## 
## loaded via a namespace (and not attached):
##  [1] httr_1.3.1                    RMySQL_0.10.14               
##  [3] bit64_0.9-7                   jsonlite_1.5                 
##  [5] AnnotationHub_2.9.17          shiny_1.0.5                  
##  [7] assertthat_0.2.0              interactiveDisplayBase_1.15.0
##  [9] blob_1.1.1                    GenomeInfoDbData_1.0.0       
## [11] Rsamtools_1.30.0              yaml_2.1.19                  
## [13] progress_1.1.2                RSQLite_2.1.0                
## [15] backports_1.1.2-9000          lattice_0.20-35              
## [17] digest_0.6.15                 XVector_0.18.0               
## [19] colorspace_1.3-2              htmltools_0.3.6              
## [21] httpuv_1.3.6.2                Matrix_1.2-10                
## [23] XML_3.98-1.10                 pkgconfig_2.0.1              
## [25] biomaRt_2.34.2                zlibbioc_1.24.0              
## [27] xtable_1.8-2                  BiocParallel_1.12.0          
## [29] SummarizedExperiment_1.8.1    lazyeval_0.2.1               
## [31] magrittr_1.5                  crayon_1.3.4                 
## [33] mime_0.5                      memoise_1.1.0                
## [35] evaluate_0.10.1               fs_1.2.2                     
## [37] MASS_7.3-47                   xml2_1.1.1                   
## [39] BiocInstaller_1.28.0          tools_3.4.1                  
## [41] prettyunits_1.0.2             matrixStats_0.52.2           
## [43] stringr_1.3.0                 DelayedArray_0.4.1           
## [45] Biostrings_2.46.0             compiler_3.4.1               
## [47] pkgdown_1.0.0                 rlang_0.2.0.9001             
## [49] grid_3.4.1                    RCurl_1.95-4.10              
## [51] htmlwidgets_1.2               bitops_1.0-6                 
## [53] rmarkdown_1.9                 curl_3.2                     
## [55] DBI_0.8                       roxygen2_6.0.1               
## [57] R6_2.2.2                      GenomicAlignments_1.14.2     
## [59] knitr_1.20                    rtracklayer_1.38.3           
## [61] bit_1.1-12                    commonmark_1.4               
## [63] rprojroot_1.3-2               ProtGenerics_1.9.1           
## [65] desc_1.2.0                    stringi_1.2.2                
## [67] Rcpp_0.12.16