The GenomicFeatures
package retrieves and manages
transcript-related features from the UCSC Genome
Bioinformatics\footnote{(http://genome.ucsc.edu/)} and
BioMart\footnote{(http://www.biomart.org/)} data resources. The
package is useful for ChIP-chip, ChIP-seq, and RNA-seq analyses.
suppressPackageStartupMessages(library('GenomicFeatures'))
TxDb
ObjectsThe GenomicFeatures
package uses TxDb
objects to store transcript metadata. This class maps the 5’ and 3’
untranslated regions (UTRs), protein coding sequences (CDSs) and exons
for a set of mRNA transcripts to their associated
genome. TxDb
objects have numerous accessors functions to
allow such features to be retrieved individually or grouped together
in a way that reflects the underlying biology.
All TxDb
objects are backed by a SQLite database that
manages genomic locations and the relationships between pre-processed
mRNA transcripts, exons, protein coding sequences, and their related
gene identifiers.
TxDb
objectsThere are two ways that users can load pre-existing data to generate a
TxDb
object. One method is to use the
loadDb
method to load the object directly from an
appropriate .sqlite database file.
Here we are loading a previously created TxDb
object
based on UCSC known gene data. This database only contains a small
subset of the possible annotations for human and is only included to
demonstrate and test the functionality of the
GenomicFeatures
package as a demonstration.
samplefile <- system.file("extdata", "hg19_knownGene_sample.sqlite",
package="GenomicFeatures")
txdb <- loadDb(samplefile)
txdb
## TxDb object:
## # Db type: TxDb
## # Supporting package: GenomicFeatures
## # Data source: UCSC
## # Genome: hg19
## # Organism: Homo sapiens
## # UCSC Table: knownGene
## # Resource URL: http://genome.ucsc.edu/
## # Type of Gene ID: Entrez Gene ID
## # Full dataset: no
## # miRBase build ID: NA
## # transcript_nrow: 178
## # exon_nrow: 620
## # cds_nrow: 523
## # Db created by: GenomicFeatures package from Bioconductor
## # Creation time: 2014-10-08 10:31:15 -0700 (Wed, 08 Oct 2014)
## # GenomicFeatures version at creation time: 1.17.21
## # RSQLite version at creation time: 0.11.4
## # DBSCHEMAVERSION: 1.0
In this case, the TxDb
object has been returned by
the loadDb
method.
More commonly however, we expect that users will just load a TxDb annotation package like this:
library(TxDb.Hsapiens.UCSC.hg19.knownGene)
txdb <- TxDb.Hsapiens.UCSC.hg19.knownGene #shorthand (for convenience)
txdb
## TxDb object:
## # Db type: TxDb
## # Supporting package: GenomicFeatures
## # Data source: UCSC
## # Genome: hg19
## # Organism: Homo sapiens
## # Taxonomy ID: 9606
## # UCSC Table: knownGene
## # Resource URL: http://genome.ucsc.edu/
## # Type of Gene ID: Entrez Gene ID
## # Full dataset: yes
## # miRBase build ID: GRCh37
## # transcript_nrow: 82960
## # exon_nrow: 289969
## # cds_nrow: 237533
## # Db created by: GenomicFeatures package from Bioconductor
## # Creation time: 2015-10-07 18:11:28 +0000 (Wed, 07 Oct 2015)
## # GenomicFeatures version at creation time: 1.21.30
## # RSQLite version at creation time: 1.0.0
## # DBSCHEMAVERSION: 1.1
Loading the package like this will also create a TxDb
object, and by default that object will have the same name as the
package itself.
It is possible to filter the data that is returned from a
TxDb
object based on it’s chromosome. This can be a
useful way to limit the things that are returned if you are only
interested in studying a handful of chromosomes.
To determine which chromosomes are currently active, use the
seqlevels
method. For example:
head(seqlevels(txdb))
## [1] "chr1" "chr2" "chr3" "chr4" "chr5" "chr6"
Will tell you all the chromosomes that are active for the
TxDb.Hsapiens.UCSC.hg19.knownGene TxDb
object (by
default it will be all of them).
If you then wanted to only set Chromosome 1 to be active you could do it like this:
seqlevels(txdb) <- "chr1"
So if you ran this, then from this point on in your R session only
chromosome 1 would be consulted when you call the various retrieval
methods… If you need to reset back to the original seqlevels (i.e.
to the seqlevels stored in the db), then set the seqlevels to
seqlevels0(txdb)
.
seqlevels(txdb) <- seqlevels0(txdb)
\begin{Exercise} Use seqlevels to set only chromsome 15 to be active. BTW, the rest of this vignette will assume you have succeeded at this. \end{Exercise}
`;N\begin{Solution}
seqlevels(txdb) <- "chr15"
seqlevels(txdb)
## [1] "chr15"
\end{Solution};N`
The TxDb
objects inherit from AnnotationDb
objects (just as the ChipDb
and OrgDb
objects do).
One of the implications of this relationship is that these object
ought to be used in similar ways to each other. Therefore we have
written supporting columns
, keytypes
, keys
and select
methods for TxDb
objects.
These methods can be a useful way of extracting data from a
TxDb
object. And they are used in the same way that
they would be used to extract information about a ChipDb
or
an OrgDb
object. Here is a simple example of how to find the
UCSC transcript names that match with a set of gene IDs.
keys <- c("100033416", "100033417", "100033420")
columns(txdb)
## [1] "CDSCHROM" "CDSEND" "CDSID" "CDSNAME" "CDSSTART"
## [6] "CDSSTRAND" "EXONCHROM" "EXONEND" "EXONID" "EXONNAME"
## [11] "EXONRANK" "EXONSTART" "EXONSTRAND" "GENEID" "TXCHROM"
## [16] "TXEND" "TXID" "TXNAME" "TXSTART" "TXSTRAND"
## [21] "TXTYPE"
keytypes(txdb)
## [1] "CDSID" "CDSNAME" "EXONID" "EXONNAME" "GENEID" "TXID" "TXNAME"
select(txdb, keys = keys, columns="TXNAME", keytype="GENEID")
## 'select()' returned 1:1 mapping between keys and columns
## GENEID TXNAME
## 1 100033416 uc001yxl.4
## 2 100033417 uc001yxo.3
## 3 100033420 uc001yxr.3
\begin{Exercise} For the genes in the example above, find the chromosome and strand information that will go with each of the transcript names. \end{Exercise} `;N\begin{Solution}
columns(txdb)
## [1] "CDSCHROM" "CDSEND" "CDSID" "CDSNAME" "CDSSTART"
## [6] "CDSSTRAND" "EXONCHROM" "EXONEND" "EXONID" "EXONNAME"
## [11] "EXONRANK" "EXONSTART" "EXONSTRAND" "GENEID" "TXCHROM"
## [16] "TXEND" "TXID" "TXNAME" "TXSTART" "TXSTRAND"
## [21] "TXTYPE"
cols <- c("TXNAME", "TXSTRAND", "TXCHROM")
select(txdb, keys=keys, columns=cols, keytype="GENEID")
## 'select()' returned 1:1 mapping between keys and columns
## GENEID TXNAME TXCHROM TXSTRAND
## 1 100033416 uc001yxl.4 chr15 +
## 2 100033417 uc001yxo.3 chr15 +
## 3 100033420 uc001yxr.3 chr15 +
\end{Solution};N`
GRanges
objectsRetrieving data with select is useful, but sometimes it is more
convenient to extract the result as GRanges
objects. This is
often the case when you are doing counting or specialized overlap
operations downstream. For these use cases there is another family of
methods available.
Perhaps the most common operations for a TxDb
object
is to retrieve the genomic coordinates or ranges for exons,
transcripts or coding sequences. The functions
transcripts
, exons
, and cds
return
the coordinate information as a GRanges
object.
As an example, all transcripts present in a TxDb
object
can be obtained as follows:
GR <- transcripts(txdb)
GR[1:3]
## GRanges object with 3 ranges and 2 metadata columns:
## seqnames ranges strand | tx_id tx_name
## <Rle> <IRanges> <Rle> | <integer> <character>
## [1] chr15 20362688-20364420 + | 53552 uc001yte.1
## [2] chr15 20487997-20496811 + | 53553 uc001ytf.1
## [3] chr15 20723929-20727150 + | 53554 uc001ytj.3
## -------
## seqinfo: 1 sequence from hg19 genome
The transcripts
function returns a GRanges
class
object. You can learn a lot more about the manipulation of these
objects by reading the GenomicRanges
introductory
vignette. The show
method for a GRanges
object
will display the ranges, seqnames (a chromosome or a contig), and
strand on the left side and then present related metadata on the right
side. At the bottom, the seqlengths display all the possible seqnames
along with the length of each sequence.
The strand
function is used to obtain the strand
information from the transcripts. The sum of the Lengths of the
Rle
object that strand
returns is equal to the
length of the GRanges
object.
tx_strand <- strand(GR)
tx_strand
## factor-Rle of length 3337 with 2 runs
## Lengths: 1732 1605
## Values : + -
## Levels(3): + - *
sum(runLength(tx_strand))
## [1] 3337
length(GR)
## [1] 3337
In addition, the transcripts
function can also be used to
retrieve a subset of the transcripts available such as those on the
\(+\)-strand of chromosome 1.
GR <- transcripts(txdb, filter=list(tx_chrom = "chr15", tx_strand = "+"))
length(GR)
## [1] 1732
unique(strand(GR))
## [1] +
## Levels: + - *
The promoters
function computes a GRanges
object
that spans the promoter region around the transcription start site
for the transcripts in a TxDb
object. The upstream
and downstream
arguments define the number of bases upstream
and downstream from the transcription start site that make up the
promoter region.
PR <- promoters(txdb, upstream=2000, downstream=400)
PR
## GRanges object with 3337 ranges and 2 metadata columns:
## seqnames ranges strand | tx_id tx_name
## <Rle> <IRanges> <Rle> | <integer> <character>
## uc001yte.1 chr15 20360688-20363087 + | 53552 uc001yte.1
## uc001ytf.1 chr15 20485997-20488396 + | 53553 uc001ytf.1
## uc001ytj.3 chr15 20721929-20724328 + | 53554 uc001ytj.3
## uc021sex.1 chr15 20737312-20739711 + | 53555 uc021sex.1
## uc010tzb.1 chr15 20740235-20742634 + | 53556 uc010tzb.1
## ... ... ... ... . ... ...
## uc021syy.1 chr15 102302656-102305055 - | 56884 uc021syy.1
## uc002cdf.1 chr15 102462863-102465262 - | 56885 uc002cdf.1
## uc002cds.2 chr15 102518897-102521296 - | 56886 uc002cds.2
## uc010utv.1 chr15 102518897-102521296 - | 56887 uc010utv.1
## uc010utw.1 chr15 102518897-102521296 - | 56888 uc010utw.1
## -------
## seqinfo: 1 sequence from hg19 genome
The exons
and cds
functions can also be used
in a similar fashion to retrive genomic coordinates for exons and
coding sequences.
\begin{Exercise} Use exonsto retrieve all the exons from chromosome 15. How does the length of this compare to the value returned by
transcripts? \end{Exercise}
`;N\begin{Solution}
EX <- exons(txdb)
EX[1:4]
## GRanges object with 4 ranges and 1 metadata column:
## seqnames ranges strand | exon_id
## <Rle> <IRanges> <Rle> | <integer>
## [1] chr15 20362688-20362858 + | 192986
## [2] chr15 20362943-20363123 + | 192987
## [3] chr15 20364397-20364420 + | 192988
## [4] chr15 20487997-20488227 + | 192989
## -------
## seqinfo: 1 sequence from hg19 genome
length(EX)
## [1] 10771
length(GR)
## [1] 1732
\end{Solution};N`
Often one is interested in how particular genomic features relate to
each other, and not just their location. For example, it might be of
interest to group transcripts by gene or to group exons by transcript.
Such groupings are supported by the transcriptsBy
,
exonsBy
, and cdsBy
functions.
The following call can be used to group transcripts by genes:
GRList <- transcriptsBy(txdb, by = "gene")
length(GRList)
## [1] 799
names(GRList)[10:13]
## [1] "100033424" "100033425" "100033427" "100033428"
GRList[11:12]
## GRangesList object of length 2:
## $`100033425`
## GRanges object with 1 range and 2 metadata columns:
## seqnames ranges strand | tx_id tx_name
## <Rle> <IRanges> <Rle> | <integer> <character>
## [1] chr15 25324204-25325381 + | 53638 uc001yxw.4
## -------
## seqinfo: 1 sequence from hg19 genome
##
## $`100033427`
## GRanges object with 1 range and 2 metadata columns:
## seqnames ranges strand | tx_id tx_name
## <Rle> <IRanges> <Rle> | <integer> <character>
## [1] chr15 25326433-25326526 + | 53640 uc001yxz.3
## -------
## seqinfo: 1 sequence from hg19 genome
The transcriptsBy
function returns a GRangesList
class object. As with GRanges
objects, you can learn more
about these objects by reading the GenomicRanges
introductory vignette. The show
method for a
GRangesList
object will display as a list of GRanges
objects. And, at the bottom the seqinfo will be displayed once for
the entire list.
For each of these three functions, there is a limited set of options
that can be passed into the by
argument to allow grouping.
For the transcriptsBy
function, you can group by gene,
exon or cds, whereas for the exonsBy
and cdsBy
functions can only be grouped by transcript (tx) or gene.
So as a further example, to extract all the exons for each transcript you can call:
GRList <- exonsBy(txdb, by = "tx")
length(GRList)
## [1] 3337
names(GRList)[10:13]
## [1] "53561" "53562" "53563" "53564"
GRList[[12]]
## GRanges object with 1 range and 3 metadata columns:
## seqnames ranges strand | exon_id exon_name exon_rank
## <Rle> <IRanges> <Rle> | <integer> <character> <integer>
## [1] chr15 22043463-22043502 + | 193028 <NA> 1
## -------
## seqinfo: 1 sequence from hg19 genome
As you can see, the GRangesList
objects returned from each
function contain locations and identifiers grouped into a list like
object according to the type of feature specified in the by
argument. The object returned can then be used by functions like
findOverlaps
to contextualize alignments from
high-throughput sequencing.
The identifiers used to label the GRanges
objects depend upon
the data source used to create the TxDb
object. So
the list identifiers will not always be Entrez Gene IDs, as they were
in the first example. Furthermore, some data sources do not provide a
unique identifier for all features. In this situation, the group
label will be a synthetic ID created by GenomicFeatures
to
keep the relations between features consistent in the database this
was the case in the 2nd example. Even though the results will
sometimes have to come back to you as synthetic IDs, you can still
always retrieve the original IDs.
\begin{Exercise} Starting with the tx_ids that are the names of the GRList object we just made, use selectto retrieve that matching transcript names. Remember that the list used a
by argument = "tx", so the list is grouped by transcript IDs. \end{Exercise}
`;N\begin{Solution}
GRList <- exonsBy(txdb, by = "tx")
tx_ids <- names(GRList)
head(select(txdb, keys=tx_ids, columns="TXNAME", keytype="TXID"))
## 'select()' returned 1:1 mapping between keys and columns
## TXID TXNAME
## 1 53552 uc001yte.1
## 2 53553 uc001ytf.1
## 3 53554 uc001ytj.3
## 4 53555 uc021sex.1
## 5 53556 uc010tzb.1
## 6 53557 uc021sey.1
\end{Solution};N`
Finally, the order of the results in a GRangesList
object can
vary with the way in which things were grouped. In most cases the
grouped elements of the GRangesList
object will be listed in
the order that they occurred along the chromosome. However, when
exons or CDS are grouped by transcript, they will instead be grouped
according to their position along the transcript itself. This is
important because alternative splicing can mean that the order along
the transcript can be different from that along the chromosome.
The intronsByTranscript
, fiveUTRsByTranscript
and threeUTRsByTranscript
are convenience functions that
provide behavior equivalent to the grouping functions, but in
prespecified form. These functions return a GRangesList
object grouped by transcript for introns, 5’ UTR’s, and 3’ UTR’s,
respectively. Below are examples of how you can call these methods.
length(intronsByTranscript(txdb))
## [1] 3337
length(fiveUTRsByTranscript(txdb))
## [1] 1825
length(threeUTRsByTranscript(txdb))
## [1] 1803
The GenomicFeatures
package also provides provides
functions for converting from ranges to actual sequence (when paired
with an appropriate BSgenome
package).
library(BSgenome.Hsapiens.UCSC.hg19)
## Loading required package: BSgenome
## Loading required package: Biostrings
## Loading required package: XVector
##
## Attaching package: 'Biostrings'
## The following object is masked from 'package:base':
##
## strsplit
## Loading required package: rtracklayer
tx_seqs1 <- extractTranscriptSeqs(Hsapiens, TxDb.Hsapiens.UCSC.hg19.knownGene,
use.names=TRUE)
And, once these sequences have been extracted, you can translate them
into proteins with translate
:
suppressWarnings(translate(tx_seqs1))
## AAStringSet object of length 3337:
## width seq names
## [1] 125 EDQDDEARVQYEGFRPGMYVRV...YTPQHMHCGAAFWA*FSDSCH uc001yte.1
## [2] 288 RIAS*GRAEFSSAQTSEIQRRR...ESVFYSVYFNYGNNCFFTVTD uc001ytf.1
## [3] 588 RSGQRLPEQPEAEGGDPGKQRR...RDLLENETHLYLCSIKICFSS uc001ytj.3
## [4] 10 HHLNCRPQTG uc021sex.1
## [5] 9 STVTLPHSQ uc010tzb.1
## ... ... ...
## [3333] 10 QVPMRVQVGQ uc021syy.1
## [3334] 306 MVTEFIFLGLSDSQELQTFLFM...DMKTAIRRLRKWDAHSSVKF* uc002cdf.1
## [3335] 550 LAVSLFFDLFFLFMCICCLLAQ...TPRRLHPAQLEILY*KHTVGF uc002cds.2
## [3336] 496 LAVSLFFDLFFLFMCICCLLAQ...PETFASCTARDPLLKAHCWFL uc010utv.1
## [3337] 531 LAVSLFFDLFFLFMCICCLLAQ...TPRRLHPAQLEILY*KHTVGF uc010utw.1
\begin{Exercise} But of course this is not a meaningful translation, because the call to extractTranscriptSeqswill have extracted all the transcribed regions of the genome regardless of whether or not they are translated. Look at the manual page for
extractTranscriptSeqs and see how you can use cdsBy to only translate only the coding regions. \end{Exercise}
`;N\begin{Solution}
cds_seqs <- extractTranscriptSeqs(Hsapiens,
cdsBy(txdb, by="tx", use.names=TRUE))
translate(cds_seqs)
## AAStringSet object of length 1875:
## width seq names
## [1] 102 MYVRVEIENVPCEFVQNIDPHY...RQRLLKYTPQHMHCGAAFWA* uc001yte.1
## [2] 435 MEWKLEQSMREQALLKAQLTQL...LGSNCCVPFFCWAWPPRRRR* uc010tzc.1
## [3] 317 MKIANNTVVTEFILLGLTQSQD...SMKRLLSRHVVCQVDFIIRN* uc001yuc.1
## [4] 314 METANYTKVTEFVLTGLSQTPE...KEVKAAMRKLVTKYILCKEK* uc010tzu.2
## [5] 317 MKIANNTVVTEFILLGLTQSQD...SMKRLLSRHVVCQVDFIIRN* uc010tzv.2
## ... ... ...
## [1871] 186 MAGGVLPLRGLRALCRVLLFLS...CLGRSEFKDICQQNVFLQVY* uc010ush.1
## [1872] 258 MYNSKLWEASGHWQHYSENMFT...PVNFLKKDLWLTLTWITVVH* uc002bxl.3
## [1873] 803 MAAEALAAEAVASRLERQEEDI...AIDKLKNLRKTRTLNAEEAF* uc002bxm.3
## [1874] 306 MVTEFIFLGLSDSQELQTFLFM...DMKTAIRRLRKWDAHSSVKF* uc002cdf.1
## [1875] 134 MSESINFSHNLGQLLSPPRCVV...KGETQESVESRVLPGPRHRH* uc010utv.1
\end{Solution};N`
TxDb
Objects or PackagesThe GenomicFeatures
package provides functions to create
TxDb
objects based on data downloaded from UCSC
Genome Bioinformatics or BioMart. The following subsections
demonstrate the use of these functions. There is also support for
creating TxDb
objects from custom data sources using
makeTxDb
; see the help page for this function for
details.
makeTxDbFromUCSC
The function makeTxDbFromUCSC
downloads UCSC
Genome Bioinformatics transcript tables (e.g. knownGene
,
refGene
, ensGene
) for a genome build (e.g.
mm9
, hg19
). Use the supportedUCSCtables
utility function to get the list of tables known to work with
makeTxDbFromUCSC
.
supportedUCSCtables(genome="mm9")
## tablename track subtrack
## 1 knownGene UCSC Genes <NA>
## 2 knownGeneOld11 Old UCSC Genes <NA>
## 3 knownGeneOld8 Old UCSC Genes <NA>
## 4 knownGeneOld7 Old UCSC Genes <NA>
## 5 knownGeneOld6 Old UCSC Genes <NA>
## 6 knownGeneOld4 Old UCSC Genes <NA>
## 7 knownGeneOld3 Old UCSC Genes <NA>
## 8 ccdsGene CCDS <NA>
## 9 refGene RefSeq Genes <NA>
## 10 xenoRefGene Other RefSeq <NA>
## 11 vegaGene Vega Genes Vega Protein Genes
## 12 vegaPseudoGene Vega Genes Vega Pseudogenes
## 13 ensGene Ensembl Genes <NA>
## 14 acembly AceView Genes <NA>
## 15 nscanPasaGene N-SCAN N-SCAN PASA-EST
## 16 nscanGene N-SCAN N-SCAN
## 17 sgpGene SGP Genes <NA>
## 18 geneid Geneid Genes <NA>
## 19 genscan Genscan Genes <NA>
## 20 exoniphy Exoniphy <NA>
## 21 augustusGene AUGUSTUS <NA>
mm9KG_txdb <- makeTxDbFromUCSC(genome="mm9", tablename="knownGene")
makeTxDbFromBiomart
Retrieve data from BioMart by specifying the mart and the data set to
the makeTxDbFromBiomart
function (not all BioMart
data sets are currently supported):
mmusculusEnsembl <- makeTxDbFromBiomart(dataset="mmusculus_gene_ensembl")
As with the makeTxDbFromUCSC
function, the
makeTxDbFromBiomart
function also has a
circ_seqs
argument that will default to using the contents
of the DEFAULT_CIRC_SEQS
vector. And just like those UCSC
sources, there is also a helper function called
getChromInfoFromBiomart
that can show what the different
chromosomes are called for a given source.
Using the makeTxDbFromBiomart
makeTxDbFromUCSC
functions can take a while and
may also require some bandwidth as these methods have to download and
then assemble a database from their respective sources. It is not
expected that most users will want to do this step every time.
Instead, we suggest that you save your annotation objects and label
them with an appropriate time stamp so as to facilitate reproducible
research.
makeTxDbFromEnsembl
The makeTxDbFromEnsembl
function creates a TxDb
object
for a given organism by importing the genomic locations of its transcripts,
exons, CDS, and genes from an Ensembl database.
See ?makeTxDbFromEnsembl
for more information.
makeTxDbFromGFF
You can also extract transcript information from either GFF3 or GTF
files by using the makeTxDbFromGFF
function.
Usage is similar to makeTxDbFromBiomart
and
makeTxDbFromUCSC
.
TxDb
ObjectOnce a TxDb
object has been created, it can be saved
to avoid the time and bandwidth costs of recreating it and to make it
possible to reproduce results with identical genomic feature data at a
later date. Since TxDb
objects are backed by a
SQLite database, the save format is a SQLite database file (which
could be accessed from programs other than R if desired). Note that
it is not possible to serialize a TxDb
object using
R’s save
function.
saveDb(mm9KG_txdb, file="fileName.sqlite")
And as was mentioned earlier, a saved TxDb
object can
be initialized from a .sqlite file by simply using loadDb
.
mm9KG_txdb <- loadDb("fileName.sqlite")
makeTxDbPackageFromUCSC
and`makeTxDbPackageFromBiomart`
It is often much more convenient to just make an annotation package
out of your annotations. If you are finding that this is the case,
then you should consider the convenience functions:
makeTxDbPackageFromUCSC
and
makeTxDbPackageFromBiomart
. These functions are similar
to makeTxDbFromUCSC
and
makeTxDbFromBiomart
except that they will take the
extra step of actually wrapping the database up into an annotation
package for you. This package can then be installed and used as of
the standard TxDb packages found on in the Bioconductor
repository.
## R version 4.3.1 (2023-06-16)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 22.04.3 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.17-bioc/R/lib/libRblas.so
## LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.10.0
##
## locale:
## [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
## [3] LC_TIME=en_GB 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
##
## time zone: America/New_York
## tzcode source: system (glibc)
##
## attached base packages:
## [1] stats4 stats graphics grDevices utils datasets methods
## [8] base
##
## other attached packages:
## [1] BSgenome.Hsapiens.UCSC.hg19_1.4.3
## [2] BSgenome_1.68.0
## [3] rtracklayer_1.60.1
## [4] Biostrings_2.68.1
## [5] XVector_0.40.0
## [6] TxDb.Hsapiens.UCSC.hg19.knownGene_3.2.2
## [7] GenomicFeatures_1.52.2
## [8] AnnotationDbi_1.62.2
## [9] Biobase_2.60.0
## [10] GenomicRanges_1.52.0
## [11] GenomeInfoDb_1.36.2
## [12] IRanges_2.34.1
## [13] S4Vectors_0.38.1
## [14] BiocGenerics_0.46.0
##
## loaded via a namespace (and not attached):
## [1] KEGGREST_1.40.0 SummarizedExperiment_1.30.2
## [3] rjson_0.2.21 xfun_0.40
## [5] lattice_0.21-8 vctrs_0.6.3
## [7] tools_4.3.1 bitops_1.0-7
## [9] generics_0.1.3 curl_5.0.2
## [11] parallel_4.3.1 tibble_3.2.1
## [13] fansi_1.0.4 RSQLite_2.3.1
## [15] blob_1.2.4 pkgconfig_2.0.3
## [17] Matrix_1.6-1 dbplyr_2.3.3
## [19] lifecycle_1.0.3 GenomeInfoDbData_1.2.10
## [21] compiler_4.3.1 stringr_1.5.0
## [23] Rsamtools_2.16.0 progress_1.2.2
## [25] codetools_0.2-19 RCurl_1.98-1.12
## [27] yaml_2.3.7 pillar_1.9.0
## [29] crayon_1.5.2 BiocParallel_1.34.2
## [31] DelayedArray_0.26.7 cachem_1.0.8
## [33] abind_1.4-5 tidyselect_1.2.0
## [35] digest_0.6.33 stringi_1.7.12
## [37] dplyr_1.1.2 restfulr_0.0.15
## [39] grid_4.3.1 biomaRt_2.56.1
## [41] fastmap_1.1.1 cli_3.6.1
## [43] magrittr_2.0.3 S4Arrays_1.0.5
## [45] XML_3.99-0.14 utf8_1.2.3
## [47] prettyunits_1.1.1 filelock_1.0.2
## [49] rappdirs_0.3.3 bit64_4.0.5
## [51] httr_1.4.7 matrixStats_1.0.0
## [53] bit_4.0.5 png_0.1-8
## [55] hms_1.1.3 memoise_2.0.1
## [57] evaluate_0.21 knitr_1.43
## [59] BiocIO_1.10.0 BiocFileCache_2.8.0
## [61] rlang_1.1.1 glue_1.6.2
## [63] DBI_1.1.3 xml2_1.3.5
## [65] R6_2.5.1 MatrixGenerics_1.12.3
## [67] GenomicAlignments_1.36.0 zlibbioc_1.46.0