This is the vignette of chimeraviz, an R package that automates the creation of chimeric RNA visualizations. This vignette will take you through the functionality of chimeraviz.
chimeraviz allows you to import data from nine different fusion-finders: deFuse, EricScript, InFusion, JAFFA, FusionCatcher, FusionMap, PRADA, SOAPFuse, and STAR-FUSION. Getting started is easy:
# Load chimeraviz
library(chimeraviz)
# Get reference to results file from deFuse
defuse833ke <- system.file(
"extdata",
"defuse_833ke_results.filtered.tsv",
package="chimeraviz")
# Load the results file into a list of fusion objects
fusions <- import_defuse(defuse833ke, "hg19")
Import functions for the other supported fusion-finders are similarly
named (for example importEriscript
or
importInfusion
).
A list of Fusion objects, objects that represent each fusion
transcript, is now available in the list fusions
.
## [1] 17
As you can see, this list has 17 fusion objects. It is straightforward to find a specific fusion event and print information about it, or about each of the partner genes.
# Find a specific fusion event
fusion <- get_fusion_by_id(fusions, 5267)
# Show information about this fusion event
fusion
## [1] "Fusion object"
## [1] "id: 5267"
## [1] "Fusion tool: defuse"
## [1] "Genome version: hg19"
## [1] "Gene names: RCC1-HENMT1"
## [1] "Chromosomes: chr1-chr1"
## [1] "Strands: +,-"
## [1] "In-frame?: FALSE"
## [1] "PartnerGene object"
## [1] "Name: RCC1"
## [1] "ensemblId: ENSG00000180198"
## [1] "Chromosome: chr1"
## [1] "Strand: +"
## [1] "Breakpoint: 28834672"
## [1] "PartnerGene object"
## [1] "Name: HENMT1"
## [1] "ensemblId: ENSG00000162639"
## [1] "Chromosome: chr1"
## [1] "Strand: -"
## [1] "Breakpoint: 109202584"
The overview plot is a nice way to get an overview over the nominated fusions in a sample. It will produce a circular plot like this one:
In this plot, you can see the following:
plot_circle()
documentation for
more detail on how this is computedAll we need for the overview plot is a list of fusion events. Let’s import the SOAPfuse example data included with chimeraviz:
# Load SOAPfuse data
if(!exists("soapfuse833ke"))
soapfuse833ke <- system.file(
"extdata",
"soapfuse_833ke_final.Fusion.specific.for.genes",
package = "chimeraviz")
fusions <- import_soapfuse(soapfuse833ke, "hg38", 10)
With these fusion events in the fusions
variable, we’re
ready to plot:
The fusion reads plot is a way to visualize the reads supporting a
fusion event mapped to the putative fusion sequence. Many fusion finders
report a putative fusion sequence, and by mapping reads to this sequence
visualize how well the fusion event is supported. The function
plot_fusion_reads()
will, given enough data, produce a plot
like this one:
As seen in the plot, this fusion event is supported by 6 paired end reads.
We have to complete a few steps in order to create the fusion reads plot. We have to:
First we load an example fusion event:
Some fusion finders report which reads that gave evidence of the
fusion event. By utilizing this information, we can avoid mapping
all reads against the fusion sequence (that can take a while).
deFuse is such a fusion finder. With the script
get_fusion_fastq.pl
that deFuse provides, we can extract
the interesting reads. For the example fusion event, we have extracted
the interesting reads into the files
reads_supporting_defuse_fusion_5267.1.fq
and
reads_supporting_defuse_fusion_5267.2.fq
. These are
included in chimeraviz
and can be referenced like this:
You can create a file containing the fusion junction sequence like this:
referenceFilename <- "reference.fa"
write_fusion_reference(fusion = fusion, filename = referenceFilename)
This will create a fasta file with the sequence.
The goal in this step is to align the reads you have against the
fusion junction sequence. There are many ways to do this.
cimeraviz
offers wrapper functions for both Bowtie and
Rsubread:
First we have to create an index for our reference sequence, the fusion junction sequence. Then we will align the reads against the reference. The steps are similar for both Bowtie and Rsubread:
# First load the bowtie functions
source(system.file(
"scripts",
"bowtie.R",
package="chimeraviz"))
# Then create index
bowtie_index(
bowtie_build_location = "/path/to/bowtie-build",
reference_fasta = reference_filename)
# And align
bowtie_align(
bowtie_location = "/path/to/bowtie",
reference_name = reference_filename,
fastq1 = fastq1,
fastq2 = fastq2,
output_bam_filename = "fusion_alignment")
The code above first loads the bowtie wrapper functions. This is
necessary because the wrapper functions are not part of the
chimeraviz
package, but are added as external scripts. We
then create an index for the reference sequence, and align the fusion
reads to the reference. The result of this will be the files
fusion_alignment.bam
,
fusion_alignment.bam.bai
, and some additional files that
Bowtie creates for indexing.
Note that you can use the parameter p
with
bowtie_align
to tell bowtie how many threads to use. The
more you have available, the faster the alignment will complete.
# First load the rsubread functions
source(system.file(
"scripts",
"rsubread.R",
package="chimeraviz"))
# Then create index
rsubread_index(reference_fasta = reference_filename)
# And align
rsubread_align(
reference_name = reference_filename,
fastq1 = fastq1,
fastq2 = fastq2,
output_bam_filename = "fusion_alignment")
The code above first loads the rsubread wrapper functions. This is
necessary because the wrapper functions are not part of the
chimeraviz
package, but are added as external scripts. We
then create an index for the reference sequence, and align the fusion
reads to the reference. The result of this will be the files
fusion_alignment.bam
,
fusion_alignment.bam.bai
, and some additional files that
Rsubread creates for indexing and aligning.
If you plan to align the original fastq files against the fusion
junction sequence, note that the alignment might take a while. The
original fastq files used to generate all chimeraviz example data were
7.3 GB each. Helping in this regard is the fact that many fusion finders
report which reads supports each fusion event. deFuse does this with the
script get_fusion_fastq.pl
(introduced in version 0.7 of
deFuse). Other fusion finders may only output which read ids supported
each fusion event. In that case, chimeraviz will help the user extract
specific read ids from the original fastq files using the function
fetch_reads_from_fastq()
. The fastq files produced by
get_fusion_fastq.pl
and
fetch_reads_from_fastq()
are much smaller, resulting in the
alignment step only taking a few seconds.
We have included an example .bam file with chimeraviz
that contains the reads supporting fusion event 5267 aligned to the
fusion junction sequence. We will continue this example with that .bam
file:
Now that we have a .bam file with the fusion reads aligned to the fusion sequence, we can import it to our fusion object:
The fusion plot is the main product of chimeraviz
,
created with the plot_fusion
function. It will create a
plot like this one:
Or, alternatively:
This plot holds a lot of information. You can see:
The fusion you can see above is the RCC1-HENMT1
fusion
described by Andreas M. Hoff et al. in the paper Identification
of Novel Fusion Genes in Testicular Germ Cell Tumors (Cancer Research,
2016).
Note that the plot reverses genes as necessary, so that the fused genes are plotted in the “correct” (5’-to-3’) relative orientation.
To build the fusion plot, we need the following:
We will go through these steps in turn.
For this example we will load data from a deFuse run. deFuse,
described in the paper deFuse:
An Algorithm for Gene Fusion Discovery in Tumor RNA-Seq Data (McPherson
et al., PloS Comput. Biol, 2011), is a Perl program that nominates
gene fusions from RNA-seq data. chimeraviz
includes example
data from a deFuse run, so let’s load it using the
import_defuse()
function:
# First find the example data
if(!exists("defuse833ke"))
defuse833ke <- system.file(
"extdata",
"defuse_833ke_results.filtered.tsv",
package = "chimeraviz")
# Then load the fusion events
fusions <- import_defuse(defuse833ke, "hg19", 1)
We should now have a list of fusion events in the
fusions
variable.
Let us try to find the RCC1-HENMT1
fusion. First, we
will search our list of fusion events for fusions involving the gene
RCC1
using the function
get_fusion_by_gene_name()
:
## [[1]]
## [1] "Fusion object"
## [1] "id: 5267"
## [1] "Fusion tool: defuse"
## [1] "Genome version: hg19"
## [1] "Gene names: RCC1-HENMT1"
## [1] "Chromosomes: chr1-chr1"
## [1] "Strands: +,-"
## [1] "In-frame?: FALSE"
As you can see, we found one fusion event involving this gene. Let us
use the get_fusion_by_id()
function to retrieve the fusion
event:
We now have the fusion event we want to plot in the variable
fusion
.
To get transcript information we will use functions from the
ensembldb
package. We will use an EnsDb
object
to retrieve transcript information for each of the partner genes in our
fusion event. Let us load an example EnsDb
object:
# First find our example EnsDb file
if(!exists("edbSqliteFile"))
edbSqliteFile <- system.file(
"extdata",
"Homo_sapiens.GRCh37.74.sqlite",
package="chimeraviz")
# Then load it
edb <- ensembldb::EnsDb(edbSqliteFile)
This example object, loaded from the file
Homo_sapiens.GRCh37.74.sqlite
that is distributed with
chimeraviz
, only holds transcript data for our specific
example fusion transcript. If you want to, you can create a full
EnsDb
by downloading Ensembl data and creating an
EnsDb
object like this:
# Create EnsDb from a downloaded .gtf file
edbSqliteFile <- ensDbFromGtf(gtf = "Homo_sapiens.GRCh37.74.gtf")
# The function above create a .sqlite file, like the one that is included with
# chimeraviz. The path to the file is stored in the edbSqliteFile variable. To
# load the database, do this:
edb <- ensembldb::EnsDb(edbSqliteFile)
If you choose to do this, note that it might take some time. You only
have to do it once, though. The .sqlite file that
ensDbFromGtf
creates can be used in later sessions, like
this:
We need a .bam file (or a .bedGraph file) if we want to include
coverage information in the fusion plot. chimeraviz
includes an example .bam file that contains only the reads in the region
of RCC1
and HENMT1
. We only need a variable
pointing to it:
With the fusion event, transcript database and bamfile ready, we can
create the plot with the function plot_fusion()
:
To collapse all transcripts into one, use the
reduce_transcripts
parameter:
plot_fusion(
fusion = fusion,
bamfile = fusion5267and11759reads,
edb = edb,
non_ucsc = TRUE,
reduce_transcripts = TRUE)
Note the non_ucsc
parameter we have set to
TRUE
. If the chromosome names in your .bam file are called
“1” instead of “chr1”, then you need to set non_ucsc = TRUE
when using plot_fusion()
.
If you have read counts in a .bedGraph file (from bedtools genomecov
for example), you can use the bedGraph
parameter instead of
the bamfile
parameter and point to a .bedGraph file.
If you are only interested in the transcripts of each partner gene in a fusion event, then it is not necessary to show all the information that is in the fusion plot. The transcripts plot is designed to give a simple view on the transcripts that might be included in the fusion transcript. It can be created like this:
# Load deFuse data
if(!exists("defuse833ke"))
defuse833ke <- system.file(
"extdata",
"defuse_833ke_results.filtered.tsv",
package = "chimeraviz")
fusions <- import_defuse(defuse833ke, "hg19", 1)
# Choose a fusion object
fusion <- get_fusion_by_id(fusions, 5267)
# Load edb
if(!exists("edbSqliteFile"))
edbSqliteFile <- system.file(
"extdata",
"Homo_sapiens.GRCh37.74.sqlite",
package="chimeraviz")
edb <- ensembldb::EnsDb(edbSqliteFile)
# Plot!
plot_transcripts(
fusion = fusion,
edb = edb)
The transcripts of the upstream gene is shown on top, and the transcripts of the downstream gene is shown on the bottom of the figure. As before, the exons that might be part of the fusion transcript are in darker colors. The interesting exons are also highlighted with grey boxes. Below each transcript, the gene name and the Ensembl transcript id is shown. This transcripts plot can also be shown with coverage data, if you are interested in seeing expression levels:
# Choose a fusion object
fusion <- get_fusion_by_id(fusions, 5267)
# Load edb
if(!exists("edbSqliteFile"))
edbSqliteFile <- system.file(
"extdata",
"Homo_sapiens.GRCh37.74.sqlite",
package="chimeraviz")
edb <- ensembldb::EnsDb(edbSqliteFile)
# Get reference to the .BAM file
if(!exists("fusion5267and11759reads"))
fusion5267and11759reads <- system.file(
"extdata",
"fusion5267and11759reads.bam",
package="chimeraviz")
# Plot!
plot_transcripts(
fusion = fusion,
edb = edb,
bamfile = fusion5267and11759reads,
non_ucsc = TRUE)
Similar to the fusion plot, it is also possible to reduce the transcripts for both genes into a single transcript track. Note the use of the parameter ylim to set a fixed limit for both y-axes in the plot:
# Load deFuse data
if(!exists("defuse833ke"))
defuse833ke <- system.file(
"extdata",
"defuse_833ke_results.filtered.tsv",
package="chimeraviz")
fusions <- import_defuse(defuse833ke, "hg19", 1)
# Choose a fusion object
fusion <- get_fusion_by_id(fusions, 5267)
# Load edb
if(!exists("edbSqliteFile"))
edbSqliteFile <- system.file(
"extdata",
"Homo_sapiens.GRCh37.74.sqlite",
package="chimeraviz")
edb <- ensembldb::EnsDb(edbSqliteFile)
# Get reference to the .BAM file
if(!exists("fusion5267and11759reads"))
fusion5267and11759reads <- system.file(
"extdata",
"fusion5267and11759reads.bam",
package="chimeraviz")
# Plot!
plot_transcripts(
fusion = fusion,
edb = edb,
bamfile = fusion5267and11759reads,
non_ucsc = TRUE,
reduce_transcripts = TRUE,
ylim = c(0,1000))
It is also possible to use a .bedGraph file instead of a .bam file in this plot.
The fusion transcript plot shows the reduced version of all exons that could be part of a fusion transcript. This is a way to view all the possible parts of a fusion transcript merged into one.
# Load deFuse data
if(!exists("defuse833ke"))
defuse833ke <- system.file(
"extdata",
"defuse_833ke_results.filtered.tsv",
package="chimeraviz")
fusions <- import_defuse(defuse833ke, "hg19", 1)
# Choose a fusion object
fusion <- get_fusion_by_id(fusions, 5267)
# Load edb
if(!exists("edbSqliteFile"))
edbSqliteFile <- system.file(
"extdata",
"Homo_sapiens.GRCh37.74.sqlite",
package="chimeraviz")
edb <- ensembldb::EnsDb(edbSqliteFile)
# Plot!
plot_fusion_transcript(
fusion,
edb)
As with previous visualizations, coverage data can be added (either with a .bam file or a .bedGraph file).
# Load deFuse data
if(!exists("defuse833ke"))
defuse833ke <- system.file(
"extdata",
"defuse_833ke_results.filtered.tsv",
package="chimeraviz")
fusions <- import_defuse(defuse833ke, "hg19", 1)
# Choose a fusion object
fusion <- get_fusion_by_id(fusions, 5267)
# Load edb
if(!exists("edbSqliteFile"))
edbSqliteFile <- system.file(
"extdata",
"Homo_sapiens.GRCh37.74.sqlite",
package="chimeraviz")
edb <- ensembldb::EnsDb(edbSqliteFile)
# Get reference to the .BAM file
if(!exists("fusion5267and11759reads"))
fusion5267and11759reads <- system.file(
"extdata",
"fusion5267and11759reads.bam",
package="chimeraviz")
# Plot!
plot_fusion_transcript(
fusion,
edb,
fusion5267and11759reads)
The fusion transcript plot merges the available transcripts of both genes, remove the parts that are not possible parts of the fusion transcript, and joins the included parts into a new transcript. With the coverage data included, it is easy to see the expression levels of the parts that might be included in a fusion transcript.
The fusion transcript plot with protein domain annotations shows a specific fusion transcript along with protein domain annotation data. If a bamfile is specified, the fusion transcript will be plotted with coverage information as well.
# Load deFuse data
if(!exists("defuse833ke"))
defuse833ke <- system.file(
"extdata",
"defuse_833ke_results.filtered.tsv",
package="chimeraviz")
fusions <- import_defuse(defuse833ke, "hg19", 1)
# Choose a fusion object
fusion <- get_fusion_by_id(fusions, 5267)
# Select transcripts
gene_upstream_transcript <- "ENST00000434290"
gene_downstream_transcript <- "ENST00000370031"
# Load edb
if(!exists("edbSqliteFile"))
edbSqliteFile <- system.file(
"extdata",
"Homo_sapiens.GRCh37.74.sqlite",
package="chimeraviz")
edb <- ensembldb::EnsDb(edbSqliteFile)
# bamfile with reads in the regions of this fusion event
if(!exists("fusion5267and11759reads"))
fusion5267and11759reads <- system.file(
"extdata",
"fusion5267and11759reads.bam",
package="chimeraviz")
# bedfile with protein domains for the transcripts in this example
if(!exists("bedfile"))
bedfile <- system.file(
"extdata",
"protein_domains_5267.bed",
package="chimeraviz")
# Plot!
plot_fusion_transcript_with_protein_domain(
fusion = fusion,
edb = edb,
bamfile = fusion5267and11759reads,
bedfile = bedfile,
gene_upstream_transcript = gene_upstream_transcript,
gene_downstream_transcript = gene_downstream_transcript,
plot_downstream_protein_domains_if_fusion_is_out_of_frame = TRUE)
Note the optional parameter
plot_downstream_protein_domains_if_fusion_is_out_of_frame
.
This let the user see protein domains encoded by the downstream gene
even though the reading frame is not kept across the fusion
junction.
This plot requires a bed file that contains the transcript level coordinates of protein domains. The bed file must be in this format:
Transcript_id | Pfam_id | Start | End | Domain_name_abbreviation | Domain_name_full |
---|---|---|---|---|---|
ENST00000434290 | PF00415 | 529 | 678 | RCC1 | Regulator of chromosome condensation (RCC1) repeat |
ENST00000434290 | PF00415 | 379 | 522 | RCC1 | Regulator of chromosome condensation (RCC1) repeat |
ENST00000434290 | PF00415 | 928 | 1056 | RCC1 | Regulator of chromosome condensation (RCC1) repeat |
ENST00000370031 | PF13489 | 289 | 753 | Methyltransf_23 | Methyltransferase domain |
All that is known about a fusion event is that a fusion-finder has scored a possible link between two genes. If there are four transcript variants of the upstream gene partner and four transcript variants of the downstream gene partner, then there are in total sixteen different splice variants of the finally processed fusion transcript. And that is only if we count the known, annotated variants of each gene. How can we make sense of all this? Plotting the transcript together as in previous plots helps, but there is a better way to visualize the putative fusion transcript: As a graph. By representing transcripts as a graph, with exons as nodes and splice junctions as edges, it is much easier to get a view on the putative fusion transcript.
# Load deFuse data
if(!exists("defuse833ke"))
defuse833ke <- system.file(
"extdata",
"defuse_833ke_results.filtered.tsv",
package="chimeraviz")
fusions <- import_defuse(defuse833ke, "hg19", 1)
# Choose a fusion object
fusion <- get_fusion_by_id(fusions, 5267)
# Plot!
plot_fusion_transcripts_graph(
fusion,
edb)
As in previous plots, exons from the upstream gene are shown in blue, and exons from the downstream gene are shown in green. This graph plot is a compact view of all possible fusion transcripts based on the known exon boundary transcripts. It is easy to see that there are four possible versions of the first exon in the upstream gene. In the downstream gene, all transcripts first use the same six exons, but they differ in the end. In total, sixteen different transcript variants for the fusion transcript are shown in the graph.
There are nine different functions that let you import fusion data, one for each fusion finder:
import_defuse()
import_ericscript()
import_fusioncatcher()
import_fusionmap()
import_infusion()
import_jaffa()
import_prada()
import_soapfuse()
import_starfusion()
Let us continue working with the deFuse example data included with chimeraviz:
# Get reference to results file from deFuse
if(!exists("defuse833ke"))
defuse833ke <- system.file(
"extdata",
"defuse_833ke_results.filtered.tsv",
package="chimeraviz")
# Load the results file into a list of fusion objects
fusions <- import_defuse(defuse833ke, "hg19")
There are three helpful functions that allow you to manage your list of fusion objects:
## [1] "Fusion object"
## [1] "id: 5267"
## [1] "Fusion tool: defuse"
## [1] "Genome version: hg19"
## [1] "Gene names: RCC1-HENMT1"
## [1] "Chromosomes: chr1-chr1"
## [1] "Strands: +,-"
## [1] "In-frame?: FALSE"
## [1] 1
## [1] 1
Creating a fusion report based on a list of fusion events is a nice
way to get an overview on the fusions in a sample. The function
create_fusion_report(fusions, "output.html")
will create an
HTML page with an overview plot and a table of searchable, sortable
fusion data. The function can be used like this:
# Load SOAPfuse data
if(!exists("soapfuse833ke"))
soapfuse833ke <- system.file(
"extdata",
"soapfuse_833ke_final.Fusion.specific.for.genes",
package = "chimeraviz")
fusions <- import_soapfuse(soapfuse833ke, "hg38", 10)
# Create report!
create_fusion_report(fusions, "output.html")
The result will be a file, output.html
, that will look
somewhat like this:
The first place to look for help is in this vignette or in the
reference manual for chimeraviz
. If the documentation and
examples there are insufficient, a good place to ask for help is on the
Bioconductor support
site. The Biostars site may
also be of help.
Any bug reports or feature ideas for chimeraviz
are very
welcome! Please report them as issues on Github.
All code in this vignette was executed in the following environment:
## R version 4.4.1 (2024-06-14)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.1 LTS
##
## Matrix products: default
## BLAS: /usr/lib/x86_64-linux-gnu/openblas-pthread/libblas.so.3
## LAPACK: /usr/lib/x86_64-linux-gnu/openblas-pthread/libopenblasp-r0.3.26.so; LAPACK version 3.12.0
##
## 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
##
## time zone: Etc/UTC
## tzcode source: system (glibc)
##
## attached base packages:
## [1] grid stats4 stats graphics grDevices utils datasets
## [8] methods base
##
## other attached packages:
## [1] chimeraviz_1.33.0 data.table_1.16.2 ensembldb_2.29.1
## [4] AnnotationFilter_1.31.0 GenomicFeatures_1.57.1 AnnotationDbi_1.69.0
## [7] Biobase_2.65.1 Gviz_1.49.0 GenomicRanges_1.57.2
## [10] Biostrings_2.73.2 GenomeInfoDb_1.41.2 XVector_0.45.0
## [13] IRanges_2.39.2 S4Vectors_0.43.2 BiocGenerics_0.51.3
## [16] BiocStyle_2.33.1
##
## loaded via a namespace (and not attached):
## [1] RColorBrewer_1.1-3 sys_3.4.3
## [3] rstudioapi_0.17.1 jsonlite_1.8.9
## [5] magrittr_2.0.3 rmarkdown_2.28
## [7] BiocIO_1.15.2 zlibbioc_1.51.2
## [9] vctrs_0.6.5 memoise_2.0.1
## [11] Rsamtools_2.21.2 RCurl_1.98-1.16
## [13] base64enc_0.1-3 htmltools_0.5.8.1
## [15] S4Arrays_1.5.11 progress_1.2.3
## [17] curl_5.2.3 SparseArray_1.5.45
## [19] Formula_1.2-5 sass_0.4.9
## [21] bslib_0.8.0 htmlwidgets_1.6.4
## [23] plyr_1.8.9 httr2_1.0.5
## [25] cachem_1.1.0 buildtools_1.0.0
## [27] GenomicAlignments_1.41.0 lifecycle_1.0.4
## [29] pkgconfig_2.0.3 Matrix_1.7-1
## [31] R6_2.5.1 fastmap_1.2.0
## [33] GenomeInfoDbData_1.2.13 MatrixGenerics_1.17.1
## [35] digest_0.6.37 colorspace_2.1-1
## [37] crosstalk_1.2.1 Hmisc_5.2-0
## [39] RSQLite_2.3.7 org.Hs.eg.db_3.20.0
## [41] filelock_1.0.3 org.Mm.eg.db_3.20.0
## [43] fansi_1.0.6 httr_1.4.7
## [45] abind_1.4-8 compiler_4.4.1
## [47] bit64_4.5.2 htmlTable_2.4.3
## [49] backports_1.5.0 BiocParallel_1.39.0
## [51] DBI_1.2.3 highr_0.11
## [53] biomaRt_2.61.3 rappdirs_0.3.3
## [55] DelayedArray_0.31.14 rjson_0.2.23
## [57] gtools_3.9.5 tools_4.4.1
## [59] foreign_0.8-87 nnet_7.3-19
## [61] glue_1.8.0 restfulr_0.0.15
## [63] checkmate_2.3.2 cluster_2.1.6
## [65] generics_0.1.3 gtable_0.3.6
## [67] BSgenome_1.73.1 hms_1.1.3
## [69] xml2_1.3.6 utf8_1.2.4
## [71] pillar_1.9.0 stringr_1.5.1
## [73] RCircos_1.2.2 dplyr_1.1.4
## [75] BiocFileCache_2.13.2 lattice_0.22-6
## [77] rtracklayer_1.65.0 bit_4.5.0
## [79] deldir_2.0-4 biovizBase_1.53.0
## [81] tidyselect_1.2.1 maketools_1.3.1
## [83] knitr_1.48 gridExtra_2.3
## [85] ProtGenerics_1.37.1 SummarizedExperiment_1.35.5
## [87] xfun_0.48 matrixStats_1.4.1
## [89] DT_0.33 stringi_1.8.4
## [91] UCSC.utils_1.1.0 lazyeval_0.2.2
## [93] yaml_2.3.10 evaluate_1.0.1
## [95] codetools_0.2-20 interp_1.1-6
## [97] tibble_3.2.1 Rgraphviz_2.49.1
## [99] graph_1.83.0 BiocManager_1.30.25
## [101] cli_3.6.3 rpart_4.1.23
## [103] munsell_0.5.1 jquerylib_0.1.4
## [105] dichromat_2.0-0.1 Rcpp_1.0.13
## [107] dbplyr_2.5.0 png_0.1-8
## [109] XML_3.99-0.17 parallel_4.4.1
## [111] ggplot2_3.5.1 blob_1.2.4
## [113] prettyunits_1.2.0 latticeExtra_0.6-30
## [115] jpeg_0.1-10 bitops_1.0-9
## [117] VariantAnnotation_1.51.2 scales_1.3.0
## [119] crayon_1.5.3 rlang_1.1.4
## [121] KEGGREST_1.45.1