This package is designed to enable the Gene Regulatory Analysis using
Variable IP (GRAVI)
workflow, as a method for detecting differential binding in ChIP-Seq
datasets. As a workflow focussed on data integration, most functions
provided by the package extraChIPs
are designed to enable
comparison across datasets. This vignette looks primarily at functions
which work with GenomicRanges
objects.
In order to use the package extraChIPs
and follow this
vignette, we recommend using the package BiocManager
hosted
on CRAN. Once this is installed, the additional packages required for
this vignette (tidyverse
, plyranges
and
Gviz
) can also be installed.
The advent of the tidyverse
has led to
tibble
objects becoming a common alternative to
data.frame
or DataFrame
objects. Simple
functions within extraChIP
enable coercion from
GRanges
, GInteractions
and
DataFrame
objects to tibble objects, with list columns
correctly handled. By default these coercion functions will coerce
GRanges
elements to a character vector. Similarly,
GRanges
represented as a character column can be coerced to
the ranges element of a GRanges
object.
First let’s coerce from a tibble
(or S3
data.frame
) to a GRanges
library(tidyverse)
library(extraChIPs)
set.seed(73)
df <- tibble(
range = c("chr1:1-10:+", "chr1:5-10:+", "chr1:5-6:+"),
gene_id = "gene1",
tx_id = paste0("transcript", 1:3),
score = runif(3)
)
df
## # A tibble: 3 × 4
## range gene_id tx_id score
## <chr> <chr> <chr> <dbl>
## 1 chr1:1-10:+ gene1 transcript1 0.442
## 2 chr1:5-10:+ gene1 transcript2 0.0831
## 3 chr1:5-6:+ gene1 transcript3 0.615
## GRanges object with 3 ranges and 3 metadata columns:
## seqnames ranges strand | gene_id tx_id score
## <Rle> <IRanges> <Rle> | <character> <character> <numeric>
## [1] chr1 1-10 + | gene1 transcript1 0.4423369
## [2] chr1 5-10 + | gene1 transcript2 0.0831099
## [3] chr1 5-6 + | gene1 transcript3 0.6146112
## -------
## seqinfo: 1 sequence from an unspecified genome; no seqlengths
Coercion back to a tibble
will place the ranges as a
character column by default. However, this can be turned off and the
conventional coercion from as.data.frame
will instead be
applied, internally wrapped in as_tibble()
## # A tibble: 3 × 4
## range gene_id tx_id score
## <chr> <chr> <chr> <dbl>
## 1 chr1:1-10:+ gene1 transcript1 0.442
## 2 chr1:5-10:+ gene1 transcript2 0.0831
## 3 chr1:5-6:+ gene1 transcript3 0.615
## # A tibble: 3 × 8
## seqnames start end width strand gene_id tx_id score
## <fct> <int> <int> <int> <fct> <chr> <chr> <dbl>
## 1 chr1 1 10 10 + gene1 transcript1 0.442
## 2 chr1 5 10 6 + gene1 transcript2 0.0831
## 3 chr1 5 6 2 + gene1 transcript3 0.615
A simple feature which may be useful for printing gene names using
rmarkdown
is contained in collapseGenes()
.
Here a character vector of gene names is collapsed into a
glue
object of length `, with gene names rendered in
italics by default.
Gene1, Gene2 and Gene3
The formation of consensus peaks incorporating ranges from multiple
replicates is a key part of many ChIP-Seq analyses. A common format
returned by tools such as mcas2 callpeak
is the
narrowPeak
(or broadPeak
) format. Sets of
these across multiple replicates can imported using the function
importPeaks()
, which returns a
GRangesList()
fl <- system.file(
c("extdata/ER_1.narrowPeak", "extdata/ER_2.narrowPeak"),
package = "extraChIPs"
)
peaks <- importPeaks(fl)
names(peaks) <- c("ER_1", "ER_2")
In the above we have loaded the peaks from two replicates from the
Estrogen Receptor (ER) in T-47D cells. To form consensus peaks we can
place a restriction on the number of replicates an overlapping peak
needs to appear in. By default, the function
makeConsensus()
sets the proportion to be zero
(p = 0)
so all peaks are retained. Note that a logical
vector/column is returned for each replicate, along with the number of
replicates the consensus peak is derived from.
## GRanges object with 14 ranges and 3 metadata columns:
## seqnames ranges strand | ER_1 ER_2 n
## <Rle> <IRanges> <Rle> | <logical> <logical> <numeric>
## [1] chr1 856458-856640 * | TRUE FALSE 1
## [2] chr1 868541-868839 * | FALSE TRUE 1
## [3] chr1 1008550-1010075 * | TRUE TRUE 2
## [4] chr1 1014770-1016015 * | TRUE TRUE 2
## [5] chr1 1051307-1051918 * | TRUE TRUE 2
## ... ... ... ... . ... ... ...
## [10] chr1 1608098-1608346 * | TRUE TRUE 2
## [11] chr1 1690460-1690641 * | FALSE TRUE 1
## [12] chr1 1790733-1790975 * | TRUE FALSE 1
## [13] chr1 1878927-1879257 * | TRUE TRUE 2
## [14] chr1 1900588-1900902 * | TRUE FALSE 1
## -------
## seqinfo: 1 sequence from an unspecified genome; no seqlengths
However, we may wish for peaks to appear in all replicates, so we can
set the argument p = 1
(or any other value 0 ≤ p ≤ 1).
## GRanges object with 6 ranges and 3 metadata columns:
## seqnames ranges strand | ER_1 ER_2 n
## <Rle> <IRanges> <Rle> | <logical> <logical> <numeric>
## [1] chr1 1008550-1010075 * | TRUE TRUE 2
## [2] chr1 1014770-1016015 * | TRUE TRUE 2
## [3] chr1 1051307-1051918 * | TRUE TRUE 2
## [4] chr1 1368372-1369009 * | TRUE TRUE 2
## [5] chr1 1608098-1608346 * | TRUE TRUE 2
## [6] chr1 1878927-1879257 * | TRUE TRUE 2
## -------
## seqinfo: 1 sequence from an unspecified genome; no seqlengths
In addition, we may wish to keep the values from the
mcols
as part of our consensus peaks, such as the
qValue
. This can be specified using the var
argument, and will return a CompressedList
as returned by
reduceMC()
, seen below.
## GRanges object with 6 ranges and 4 metadata columns:
## seqnames ranges strand | qValue ER_1 ER_2
## <Rle> <IRanges> <Rle> | <NumericList> <logical> <logical>
## [1] chr1 1008550-1010075 * | 635.748,605.040 TRUE TRUE
## [2] chr1 1014770-1016015 * | 95.1223,78.3953 TRUE TRUE
## [3] chr1 1051307-1051918 * | 65.5538,95.9677 TRUE TRUE
## [4] chr1 1368372-1369009 * | 22.5483,29.1231 TRUE TRUE
## [5] chr1 1608098-1608346 * | 72.3235,64.8439 TRUE TRUE
## [6] chr1 1878927-1879257 * | 33.8900,14.4909 TRUE TRUE
## n
## <numeric>
## [1] 2
## [2] 2
## [3] 2
## [4] 2
## [5] 2
## [6] 2
## -------
## seqinfo: 1 sequence from an unspecified genome; no seqlengths
We could even identify the peak centre and pass these to the set of consensus peaks for downstream peak re-centreing.
library(plyranges)
peaks %>%
endoapply(mutate, centre = start + peak) %>%
makeConsensus(p = 1, var = "centre")
## GRanges object with 6 ranges and 4 metadata columns:
## seqnames ranges strand | centre ER_1 ER_2
## <Rle> <IRanges> <Rle> | <NumericList> <logical> <logical>
## [1] chr1 1008550-1010075 * | 1009212,1009212 TRUE TRUE
## [2] chr1 1014770-1016015 * | 1014981,1014949 TRUE TRUE
## [3] chr1 1051307-1051918 * | 1051565,1051634 TRUE TRUE
## [4] chr1 1368372-1369009 * | 1368613,1368767 TRUE TRUE
## [5] chr1 1608098-1608346 * | 1608247,1608247 TRUE TRUE
## [6] chr1 1878927-1879257 * | 1879087,1879094 TRUE TRUE
## n
## <numeric>
## [1] 2
## [2] 2
## [3] 2
## [4] 2
## [5] 2
## [6] 2
## -------
## seqinfo: 1 sequence from an unspecified genome; no seqlengths
We could then find the mean of the peak centres for an averaged centre.
peaks %>%
endoapply(mutate, centre = start + peak) %>%
makeConsensus(p = 1, var = "centre") %>%
mutate(centre = vapply(centre, mean, numeric(1)))
## GRanges object with 6 ranges and 4 metadata columns:
## seqnames ranges strand | centre ER_1 ER_2 n
## <Rle> <IRanges> <Rle> | <numeric> <logical> <logical> <numeric>
## [1] chr1 1008550-1010075 * | 1009212 TRUE TRUE 2
## [2] chr1 1014770-1016015 * | 1014965 TRUE TRUE 2
## [3] chr1 1051307-1051918 * | 1051600 TRUE TRUE 2
## [4] chr1 1368372-1369009 * | 1368690 TRUE TRUE 2
## [5] chr1 1608098-1608346 * | 1608247 TRUE TRUE 2
## [6] chr1 1878927-1879257 * | 1879090 TRUE TRUE 2
## -------
## seqinfo: 1 sequence from an unspecified genome; no seqlengths
mcols()
The standard set operations implemented in the package
GenomicRanges
will always drop the mcols
element by default. The extraChIPs
functions
reduceMC()
, setdiffMC()
,
intersectMC()
and unionMC()
all produce the
same output as their similarly-named functions, however, the
mcols()
elements in the query object are also returned.
Where required, columns are coerced into CompressedList
columns. This can particularly useful when needed to propagate the
information contained in the initial ranges through to subsequent
analytic steps
GRanges
objectsThe classical approach to defining TSS regions for a set of
transcripts would be to use the function resize
()`, setting
the width as 1.
## GRanges object with 3 ranges and 3 metadata columns:
## seqnames ranges strand | gene_id tx_id score
## <Rle> <IRanges> <Rle> | <character> <character> <numeric>
## [1] chr1 1 + | gene1 transcript1 0.4423369
## [2] chr1 5 + | gene1 transcript2 0.0831099
## [3] chr1 5 + | gene1 transcript3 0.6146112
## -------
## seqinfo: 1 sequence from an unspecified genome; no seqlengths
As we can see, two transcripts start at position 5, so we may choose
to reduce this, which would lose the mcols
element. The
alternative reduceMC()
will retain all
mcols
.
## GRanges object with 2 ranges and 0 metadata columns:
## seqnames ranges strand
## <Rle> <IRanges> <Rle>
## [1] chr1 1 +
## [2] chr1 5 +
## -------
## seqinfo: 1 sequence from an unspecified genome; no seqlengths
## GRanges object with 2 ranges and 3 metadata columns:
## seqnames ranges strand | gene_id tx_id
## <Rle> <IRanges> <Rle> | <character> <CharacterList>
## [1] chr1 1 + | gene1 transcript1
## [2] chr1 5 + | gene1 transcript2,transcript3
## score
## <NumericList>
## [1] 0.442337
## [2] 0.0831099,0.6146112
## -------
## seqinfo: 1 sequence from an unspecified genome; no seqlengths
By default, this function will attempt to coerce mcols
to a new mcol
of the appropriate type, however, when
multiple values are inevitable such as for the tx_id
column
above, these will be coerced to a CompressedList
. The
simplification of the multiple values seen in the gene_id
can also be turned off if desired should repeated values be important
for downstream analysis.
## GRanges object with 2 ranges and 3 metadata columns:
## seqnames ranges strand | gene_id tx_id
## <Rle> <IRanges> <Rle> | <CharacterList> <CharacterList>
## [1] chr1 1 + | gene1 transcript1
## [2] chr1 5 + | gene1,gene1 transcript2,transcript3
## score
## <NumericList>
## [1] 0.442337
## [2] 0.0831099,0.6146112
## -------
## seqinfo: 1 sequence from an unspecified genome; no seqlengths
This allows for simple integration with tidyverse
nesting strategies.
## # A tibble: 3 × 4
## range gene_id tx_id score
## <chr> <chr> <chr> <dbl>
## 1 chr1:1:+ gene1 transcript1 0.442
## 2 chr1:5:+ gene1 transcript2 0.0831
## 3 chr1:5:+ gene1 transcript3 0.615
Whilst reduceMC
relies on the range-reduction as
implemented in GenomicRanges::reduce()
, some alternative
approaches are included, such as chopMC()
, which finds
identical ranges and nests the mcols element as
CompressedList
objects.
## GRanges object with 2 ranges and 3 metadata columns:
## seqnames ranges strand | gene_id tx_id
## <Rle> <IRanges> <Rle> | <character> <CharacterList>
## [1] chr1 1 + | gene1 transcript1
## [2] chr1 5 + | gene1 transcript2,transcript3
## score
## <NumericList>
## [1] 0.442337
## [2] 0.0831099,0.6146112
## -------
## seqinfo: 1 sequence from an unspecified genome; no seqlengths
In the case of the object tss
, this output is identical
to reduceMC()
, however, given that there are no identical
ranges in gr
the two functions would behave very
differently on that object.
A final operation for simplifying GRanges
objects would
be distinctMC()
which is a wrapper to
dplyr]::distinct
incorporating both the range and
mcols
. The columns to search can be called using
<data-masking>
approaches as detailed in the
manual.
## GRanges object with 2 ranges and 0 metadata columns:
## seqnames ranges strand
## <Rle> <IRanges> <Rle>
## [1] chr1 1 +
## [2] chr1 5 +
## -------
## seqinfo: 1 sequence from an unspecified genome; no seqlengths
## GRanges object with 2 ranges and 1 metadata column:
## seqnames ranges strand | gene_id
## <Rle> <IRanges> <Rle> | <character>
## [1] chr1 1 + | gene1
## [2] chr1 5 + | gene1
## -------
## seqinfo: 1 sequence from an unspecified genome; no seqlengths
GRanges
ObjectsWhilst reduce/reduceMC
is applied to a single
GRanges
object, the set operation functions
intersect
, setdiff
and union
are
valuable approaches for comparing ranges. Using the *MC()
functions will retain mcols
elements from the query
range.
peaks <- GRanges(c("chr1:1-5", "chr1:9-12:*"))
peaks$peak_id <- c("peak1", "peak2")
GenomicRanges::intersect(gr, peaks, ignore.strand = TRUE)
## GRanges object with 2 ranges and 0 metadata columns:
## seqnames ranges strand
## <Rle> <IRanges> <Rle>
## [1] chr1 1-5 *
## [2] chr1 9-10 *
## -------
## seqinfo: 1 sequence from an unspecified genome; no seqlengths
## GRanges object with 2 ranges and 3 metadata columns:
## seqnames ranges strand | gene_id
## <Rle> <IRanges> <Rle> | <character>
## [1] chr1 1-5 * | gene1
## [2] chr1 9-10 * | gene1
## tx_id score
## <CharacterList> <NumericList>
## [1] transcript1,transcript2,transcript3 0.4423369,0.0831099,0.6146112
## [2] transcript1,transcript2 0.4423369,0.0831099
## -------
## seqinfo: 1 sequence from an unspecified genome; no seqlengths
## GRanges object with 1 range and 3 metadata columns:
## seqnames ranges strand | gene_id
## <Rle> <IRanges> <Rle> | <character>
## [1] chr1 6-8 * | gene1
## tx_id score
## <CharacterList> <NumericList>
## [1] transcript1,transcript2,transcript3 0.4423369,0.0831099,0.6146112
## -------
## seqinfo: 1 sequence from an unspecified genome; no seqlengths
## GRanges object with 1 range and 3 metadata columns:
## seqnames ranges strand | gene_id
## <Rle> <IRanges> <Rle> | <character>
## [1] chr1 1-12 * | gene1
## tx_id score
## <CharacterList> <NumericList>
## [1] transcript1,transcript2,transcript3 0.4423369,0.0831099,0.6146112
## -------
## seqinfo: 1 sequence from an unspecified genome; no seqlengths
There is a performance overhead to preparation of mcols as
CompressedList
objects and all mcols
in the
query object will be returned. If only wishing to retain a subset of
mcols
, these should be selected prior to passing to these
functions.
## GRanges object with 2 ranges and 1 metadata column:
## seqnames ranges strand | tx_id
## <Rle> <IRanges> <Rle> | <CharacterList>
## [1] chr1 1-5 * | transcript1,transcript2,transcript3
## [2] chr1 9-10 * | transcript1,transcript2
## -------
## seqinfo: 1 sequence from an unspecified genome; no seqlengths
Whilst the functions findOverlaps()
and
overlapsAny()
are extremely useful, the addition of
propOverlap()
returns a numeric vector with the proportion
of each query range (x
) which overlaps any range in the
subject range (y
)
## [1] 0.7 0.5 0.5
This is also extended to enable comparison across multiple features and to classify each peak by which features that it overlaps the most strongly.
## [1] "peak1" "peak2" "peak1"
In addition to standard set operations, one set of ranges can be used
to partition another set of ranges, returning mcols
from
both ranges. Ranges from the query range (x
) are returned
after being partitioned by the ranges in the subject range
(y
). Subject ranges used for partitioning must be
non-overlapping, and if overlapping ranges are provided, these will
be reduced prior to partitioning.
This enables the identification of specific ranges from the query
range (x
) which overlap ranges from the subject range
(y
) Under this approach, mcols
from both
query
and subject
ranges will be returned to
enable the clear ranges which are common and distinct within the two
sets of ranges.
## GRanges object with 5 ranges and 4 metadata columns:
## seqnames ranges strand | peak_id gene_id
## <Rle> <IRanges> <Rle> | <character> <character>
## [1] chr1 1-5 + | peak1 gene1
## [2] chr1 5 + | peak1 gene1
## [3] chr1 6 + | <NA> gene1
## [4] chr1 6-8 + | <NA> gene1
## [5] chr1 9-10 + | peak2 gene1
## tx_id score
## <CharacterList> <NumericList>
## [1] transcript1 0.442337
## [2] transcript2,transcript3 0.0831099,0.6146112
## [3] transcript3 0.614611
## [4] transcript1,transcript2 0.4423369,0.0831099
## [5] transcript1,transcript2 0.4423369,0.0831099
## -------
## seqinfo: 1 sequence from an unspecified genome; no seqlengths
Whilst this shares some similarity with intersectMC()
the additional capabilities provide greater flexibility.
## GRanges object with 2 ranges and 4 metadata columns:
## seqnames ranges strand | peak_id gene_id
## <Rle> <IRanges> <Rle> | <character> <character>
## [1] chr1 6 + | <NA> gene1
## [2] chr1 6-8 + | <NA> gene1
## tx_id score
## <CharacterList> <NumericList>
## [1] transcript3 0.614611
## [2] transcript1,transcript2 0.4423369,0.0831099
## -------
## seqinfo: 1 sequence from an unspecified genome; no seqlengths
Using the function stitchRanges()
we are able to join
together sets of nearby ranges, but with the option of placing clear
barriers between ranges, across which ranges cannot be joined. This may
be useful if joining enhancers to form putative super-enhancers, but
explicitly omitting defined promoter regions.
enh <- GRanges(c("chr1:1-10", "chr1:101-110", "chr1:181-200"))
prom <- GRanges("chr1:150:+")
se <- stitchRanges(enh, exclude = prom, maxgap = 100)
se
## GRanges object with 2 ranges and 0 metadata columns:
## seqnames ranges strand
## <Rle> <IRanges> <Rle>
## [1] chr1 1-110 *
## [2] chr1 181-200 *
## -------
## seqinfo: 1 sequence from an unspecified genome; no seqlengths
As a visualisation (below) ranges within 100bp were stitched together, however the region defined as a ‘promoter’ acted as a barrier and ranges were not stitched together across this barrier.
## 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] stats4 stats graphics grDevices utils datasets methods
## [8] base
##
## other attached packages:
## [1] magrittr_2.0.3 quantro_1.39.0
## [3] here_1.0.1 ggrepel_0.9.6
## [5] glue_1.8.0 scales_1.3.0
## [7] plyranges_1.25.0 extraChIPs_1.11.0
## [9] ggside_0.3.1 patchwork_1.3.0
## [11] edgeR_4.3.21 limma_3.61.12
## [13] rtracklayer_1.65.0 BiocParallel_1.41.0
## [15] csaw_1.39.0 SummarizedExperiment_1.35.5
## [17] Biobase_2.67.0 MatrixGenerics_1.17.1
## [19] matrixStats_1.4.1 Rsamtools_2.21.2
## [21] Biostrings_2.75.0 XVector_0.45.0
## [23] GenomicRanges_1.57.2 GenomeInfoDb_1.41.2
## [25] IRanges_2.39.2 S4Vectors_0.43.2
## [27] BiocGenerics_0.53.0 lubridate_1.9.3
## [29] forcats_1.0.0 stringr_1.5.1
## [31] dplyr_1.1.4 purrr_1.0.2
## [33] readr_2.1.5 tidyr_1.3.1
## [35] tibble_3.2.1 ggplot2_3.5.1
## [37] tidyverse_2.0.0 BiocStyle_2.35.0
##
## loaded via a namespace (and not attached):
## [1] ProtGenerics_1.37.1 bitops_1.0-9
## [3] httr_1.4.7 RColorBrewer_1.1-3
## [5] doParallel_1.0.17 InteractionSet_1.33.0
## [7] tools_4.4.1 doRNG_1.8.6
## [9] backports_1.5.0 utf8_1.2.4
## [11] R6_2.5.1 HDF5Array_1.33.8
## [13] mgcv_1.9-1 lazyeval_0.2.2
## [15] Gviz_1.49.0 rhdf5filters_1.17.0
## [17] withr_3.0.2 prettyunits_1.2.0
## [19] gridExtra_2.3 base64_2.0.2
## [21] VennDiagram_1.7.3 preprocessCore_1.67.1
## [23] cli_3.6.3 formatR_1.14
## [25] labeling_0.4.3 sass_0.4.9
## [27] genefilter_1.87.0 askpass_1.2.1
## [29] foreign_0.8-87 siggenes_1.79.0
## [31] illuminaio_0.47.0 rentrez_1.2.3
## [33] dichromat_2.0-0.1 scrime_1.3.5
## [35] BSgenome_1.75.0 rstudioapi_0.17.1
## [37] RSQLite_2.3.7 generics_0.1.3
## [39] BiocIO_1.17.0 Matrix_1.7-1
## [41] interp_1.1-6 futile.logger_1.4.3
## [43] fansi_1.0.6 abind_1.4-8
## [45] lifecycle_1.0.4 yaml_2.3.10
## [47] rhdf5_2.49.0 SparseArray_1.5.45
## [49] BiocFileCache_2.15.0 grid_4.4.1
## [51] blob_1.2.4 crayon_1.5.3
## [53] lattice_0.22-6 ComplexUpset_1.3.3
## [55] GenomicFeatures_1.57.1 annotate_1.85.0
## [57] KEGGREST_1.45.1 sys_3.4.3
## [59] maketools_1.3.1 pillar_1.9.0
## [61] knitr_1.48 beanplot_1.3.1
## [63] metapod_1.13.0 rjson_0.2.23
## [65] codetools_0.2-20 data.table_1.16.2
## [67] vctrs_0.6.5 png_0.1-8
## [69] gtable_0.3.6 cachem_1.1.0
## [71] xfun_0.48 S4Arrays_1.5.11
## [73] survival_3.7-0 iterators_1.0.14
## [75] statmod_1.5.0 GenomicInteractions_1.39.0
## [77] nlme_3.1-166 bit64_4.5.2
## [79] progress_1.2.3 filelock_1.0.3
## [81] rprojroot_2.0.4 bslib_0.8.0
## [83] nor1mix_1.3-3 rpart_4.1.23
## [85] colorspace_2.1-1 DBI_1.2.3
## [87] Hmisc_5.2-0 nnet_7.3-19
## [89] tidyselect_1.2.1 bit_4.5.0
## [91] compiler_4.4.1 curl_5.2.3
## [93] httr2_1.0.5 htmlTable_2.4.3
## [95] xml2_1.3.6 DelayedArray_0.33.1
## [97] checkmate_2.3.2 quadprog_1.5-8
## [99] rappdirs_0.3.3 digest_0.6.37
## [101] rmarkdown_2.28 GEOquery_2.73.5
## [103] htmltools_0.5.8.1 pkgconfig_2.0.3
## [105] jpeg_0.1-10 base64enc_0.1-3
## [107] sparseMatrixStats_1.17.2 highr_0.11
## [109] dbplyr_2.5.0 fastmap_1.2.0
## [111] ensembldb_2.29.1 rlang_1.1.4
## [113] htmlwidgets_1.6.4 UCSC.utils_1.1.0
## [115] DelayedMatrixStats_1.29.0 farver_2.1.2
## [117] jquerylib_0.1.4 jsonlite_1.8.9
## [119] mclust_6.1.1 VariantAnnotation_1.51.2
## [121] RCurl_1.98-1.16 Formula_1.2-5
## [123] GenomeInfoDbData_1.2.13 Rhdf5lib_1.27.0
## [125] munsell_0.5.1 Rcpp_1.0.13
## [127] stringi_1.8.4 zlibbioc_1.51.2
## [129] MASS_7.3-61 plyr_1.8.9
## [131] bumphunter_1.49.0 minfi_1.51.0
## [133] parallel_4.4.1 deldir_2.0-4
## [135] splines_4.4.1 multtest_2.61.0
## [137] hms_1.1.3 locfit_1.5-9.10
## [139] igraph_2.1.1 rngtools_1.5.2
## [141] buildtools_1.0.0 biomaRt_2.63.0
## [143] futile.options_1.0.1 XML_3.99-0.17
## [145] evaluate_1.0.1 latticeExtra_0.6-30
## [147] biovizBase_1.55.0 lambda.r_1.2.4
## [149] BiocManager_1.30.25 tzdb_0.4.0
## [151] foreach_1.5.2 tweenr_2.0.3
## [153] openssl_2.2.2 polyclip_1.10-7
## [155] reshape_0.8.9 ggforce_0.4.2
## [157] broom_1.0.7 xtable_1.8-4
## [159] restfulr_0.0.15 AnnotationFilter_1.31.0
## [161] memoise_2.0.1 AnnotationDbi_1.69.0
## [163] GenomicAlignments_1.41.0 cluster_2.1.6
## [165] timechange_0.3.0