library(ggplot2)
library(GenomicRanges)
## Loading required package: stats4
## Loading required package: BiocGenerics
## Loading required package: generics
##
## Attaching package: 'generics'
## The following objects are masked from 'package:base':
##
## as.difftime, as.factor, as.ordered, intersect, is.element,
## setdiff, setequal, union
##
## Attaching package: 'BiocGenerics'
## The following objects are masked from 'package:generics':
##
## intersect, setdiff, setequal, union
## The following objects are masked from 'package:stats':
##
## IQR, mad, sd, var, xtabs
## The following objects are masked from 'package:base':
##
## Filter, Find, Map, Position, Reduce, anyDuplicated, aperm,
## append, as.data.frame, basename, cbind, colnames, dirname,
## do.call, duplicated, eval, evalq, get, grep, grepl, intersect,
## is.unsorted, lapply, mapply, match, mget, order, paste, pmax,
## pmax.int, pmin, pmin.int, rank, rbind, rownames, sapply,
## saveRDS, setdiff, setequal, table, tapply, union, unique,
## unsplit, which.max, which.min
## Loading required package: S4Vectors
##
## Attaching package: 'S4Vectors'
## The following object is masked from 'package:utils':
##
## findMatches
## The following objects are masked from 'package:base':
##
## I, expand.grid, unname
## Loading required package: IRanges
## Loading required package: GenomeInfoDb
library(InteractionSet)
## Loading required package: SummarizedExperiment
## Loading required package: MatrixGenerics
## Loading required package: matrixStats
##
## Attaching package: 'MatrixGenerics'
## The following objects are masked from 'package:matrixStats':
##
## colAlls, colAnyNAs, colAnys, colAvgsPerRowSet, colCollapse,
## colCounts, colCummaxs, colCummins, colCumprods, colCumsums,
## colDiffs, colIQRDiffs, colIQRs, colLogSumExps, colMadDiffs,
## colMads, colMaxs, colMeans2, colMedians, colMins, colOrderStats,
## colProds, colQuantiles, colRanges, colRanks, colSdDiffs, colSds,
## colSums2, colTabulates, colVarDiffs, colVars, colWeightedMads,
## colWeightedMeans, colWeightedMedians, colWeightedSds,
## colWeightedVars, rowAlls, rowAnyNAs, rowAnys, rowAvgsPerColSet,
## rowCollapse, rowCounts, rowCummaxs, rowCummins, rowCumprods,
## rowCumsums, rowDiffs, rowIQRDiffs, rowIQRs, rowLogSumExps,
## rowMadDiffs, rowMads, rowMaxs, rowMeans2, rowMedians, rowMins,
## rowOrderStats, rowProds, rowQuantiles, rowRanges, rowRanks,
## rowSdDiffs, rowSds, rowSums2, rowTabulates, rowVarDiffs,
## rowVars, rowWeightedMads, rowWeightedMeans, rowWeightedMedians,
## rowWeightedSds, rowWeightedVars
## Loading required package: Biobase
## Welcome to Bioconductor
##
## Vignettes contain introductory material; view with
## 'browseVignettes()'. To cite Bioconductor, see
## 'citation("Biobase")', and for packages 'citation("pkgname")'.
##
## Attaching package: 'Biobase'
## The following object is masked from 'package:MatrixGenerics':
##
## rowMedians
## The following objects are masked from 'package:matrixStats':
##
## anyMissing, rowMedians
library(HiCExperiment)
## Consider using the `HiContacts` package to perform advanced genomic operations
## on `HiCExperiment` objects.
##
## Read "Orchestrating Hi-C analysis with Bioconductor" online book to learn more:
## https://js2264.github.io/OHCA/
##
## Attaching package: 'HiCExperiment'
## The following object is masked from 'package:SummarizedExperiment':
##
## metadata<-
## The following object is masked from 'package:S4Vectors':
##
## metadata<-
## The following object is masked from 'package:ggplot2':
##
## resolution
library(HiContactsData)
## Loading required package: ExperimentHub
## Loading required package: AnnotationHub
## Loading required package: BiocFileCache
## Loading required package: dbplyr
##
## Attaching package: 'AnnotationHub'
## The following object is masked from 'package:Biobase':
##
## cache
library(HiContacts)
## Registered S3 methods overwritten by 'readr':
## method from
## as.data.frame.spec_tbl_df vroom
## as_tibble.spec_tbl_df vroom
## format.col_spec vroom
## print.col_spec vroom
## print.collector vroom
## print.date_names vroom
## print.locale vroom
## str.col_spec vroom
library(rtracklayer)
##
## Attaching package: 'rtracklayer'
## The following object is masked from 'package:AnnotationHub':
##
## hubUrl
5 Matrix-centric analysis
This chapter focuses on the various analytical tools offered by HiContacts
to compute matrix-related metrics from a HiCExperiment
object.
In the first part of this book, we have seen how to query parts or all of the data contained in Hi-C contact matrices using the HiCExperiment
object (Chapter 2), how to manipulate HiCExperiment
objects (Chapter 3) and how to visualize Hi-C contact matrices as heatmaps (Chapter 4).
The HiContacts
package directly operates on HiCExperiment
objects and extends its usability by providing a comprehensive toolkit to analyze Hi-C data, focusing on four main topics:
- Contact matrix-centric analyses (this chapter)
- Interactions-centric analyses (Chapter 6)
- Structural feature annotations (Chapter 7)
- Hi-C visualization (see previous chapter)
Matrix-centric analyses consider a HiCExperiment
object from the βmatrixβ perspective to perform a range of matrix-based operations. This encompasses:
- Computing observed/expected (O/E) map
- Computing auto-correlation map
- Smoothing out a contact map
- Merging multiple Hi-C maps together
- Comparing two Hi-C maps to each other
- All the functions described in this chapter are endomorphisms: they take
HiCExperiment
objects as input and return modifiedHiCExperiment
objects. - Internally, most of the functions presented in this chapter make a call to
as.matrix(<HiCExperiment>)
to coerce it into amatrix
.
hic
object π
To demonstrate HiContacts
functionalities, we will create an HiCExperiment
object from an example .cool
file provided in the HiContactsData
package.
library(HiCExperiment)
library(HiContactsData)
# ---- This downloads an example `.mcool` file and caches it locally
coolf <- HiContactsData('yeast_wt', 'mcool')
## see ?HiContactsData and browseVignettes('HiContactsData') for documentation
## loading from cache
# ---- This creates a connection to the disk-stored `.mcool` file
cf <- CoolFile(coolf)
cf
## CoolFile object
## .mcool file: /home/biocbuild/.cache/R/ExperimentHub/1a8be61e27f712_7752
## resolution: 1000
## pairs file:
## metadata(0):
# ---- This imports contacts from the chromosome `II` at resolution `2000`
hic <- import(cf, focus = 'II', resolution = 2000)
hic
## `HiCExperiment` object with 471,364 contacts over 407 regions
## -------
## fileName: "/home/biocbuild/.cache/R/ExperimentHub/1a8be61e27f712_7752"
## focus: "II"
## resolutions(5): 1000 2000 4000 8000 16000
## active resolution: 2000
## interactions: 34063
## scores(2): count balanced
## topologicalFeatures: compartments(0) borders(0) loops(0) viewpoints(0)
## pairsFile: N/A
## metadata(0):
5.1 Operations in an individual matrix
5.1.1 Balancing a raw interaction count map
Hi-C sequencing coverage is systematically affected by multiple confounding factors, e.g. density of restriction sites, GC%, genome mappability, etc.. Overall, it generally ends up not homogenous throughout the entire genome and this leads to artifacts in un-normalized count
matrices.
To correct for sequencing coverage heterogeneity of raw count
maps, Hi-C data can be normalized using matrix balancing approaches (Cournac et al. (2012), Imakaev et al. (2012)). This is generally done directly on the disk-stored matrices using out-of-memory strategies (e.g. with cooler balance <.cool>
). However, if contact matrix files are imported into a HiCExperiment
object but no balanced
scores are available, in-memory balancing can be performed using the normalize
function. This adds an extra ICE
element in scores
list (while the interactions
themselves are unmodified).
normalized_hic <- normalize(hic)
normalized_hic
## `HiCExperiment` object with 471,364 contacts over 407 regions
## -------
## fileName: "/home/biocbuild/.cache/R/ExperimentHub/1a8be61e27f712_7752"
## focus: "II"
## resolutions(5): 1000 2000 4000 8000 16000
## active resolution: 2000
## interactions: 34063
## scores(3): count balanced ICE
## topologicalFeatures: compartments(0) borders(0) loops(0) viewpoints(0)
## pairsFile: N/A
## metadata(0):
It is possible to plot the different scores
of the resulting object to visualize the newly computed scores
. In this example, ICE
scores should be nearly identical to balanced
scores, which were originally imported from the disk-stored contact matrix.
cowplot::plot_grid(
plotMatrix(normalized_hic, use.scores = 'count', caption = FALSE),
plotMatrix(normalized_hic, use.scores = 'balanced', caption = FALSE),
plotMatrix(normalized_hic, use.scores = 'ICE', caption = FALSE),
nrow = 1
)
5.1.2 Computing observed/expected (O/E) map
The most prominent feature of a balanced Hi-C matrix is the strong main diagonal. This main diagonal is observed because interactions between immediate adjacent genomic loci are more prone to happen than interactions spanning longer genomic distances. This βexpectedβ behavior is due to the polymer nature of the chromosomes being studied, and can be locally estimated using the distance-dependent interaction frequency (a.k.a. the βdistance lawβ, or P(s)). It can be used to compute an expected
matrix on interactions.
When it is desirable to βmaskβ this polymer behavior to emphasize topological structures formed by chromosomes, one can divide a given balanced matrix by its expected
matrix, i.e. calculate the observed/expected (O/E) map. This is sometimes called βdetrendingβ, as it effectively removes the average polymer behavior from the balanced matrix.
The detrend
function performs this operation on a given HiCExperiment
object. It adds two extra elements in scores
list: expected
and detrended
metrics (while the interactions
themselves are unmodified).
detrended_hic <- detrend(hic)
detrended_hic
## `HiCExperiment` object with 471,364 contacts over 407 regions
## -------
## fileName: "/home/biocbuild/.cache/R/ExperimentHub/1a8be61e27f712_7752"
## focus: "II"
## resolutions(5): 1000 2000 4000 8000 16000
## active resolution: 2000
## interactions: 34063
## scores(4): count balanced expected detrended
## topologicalFeatures: compartments(0) borders(0) loops(0) viewpoints(0)
## pairsFile: N/A
## metadata(0):
Topological features will be visually more prominent in the O/E detrended
Hi-C map.
cowplot::plot_grid(
plotMatrix(detrended_hic, use.scores = 'balanced', scale = 'log10', limits = c(-3.5, -1.2), caption = FALSE),
plotMatrix(detrended_hic, use.scores = 'expected', scale = 'log10', limits = c(-3.5, -1.2), caption = FALSE),
plotMatrix(detrended_hic, use.scores = 'detrended', scale = 'linear', limits = c(-1, 1), cmap = bwrColors(), caption = FALSE),
nrow = 1
)
detrended
scores
-
expected
scores are inlinear
scale and Β± in the same amplitude thanbalanced
scores; -
detrended
scores are inlog2
scale, in general approximately centered around 0. When plottingdetrended
scores,scale = linear
should be set to prevent the defaultlog10
scaling.
5.1.4 Despeckling (smoothing out) a contact map
Shallow-sequenced Hi-C libraries or matrices binned with an overly small bin size sometimes produce βgrainyβ Hi-C maps with noisy backgrounds. A grainy map may also be obtained when dividing two matrices, e.g. when computing the O/E ratio with detrend
. This is particularly true for sparser long-range interactions. To overcome such limitations, HiCExperiment
objects can be βdespeckle
dβ to smooth out focal speckles.
hic2 <- detrend(hic['II:400000-700000'])
hic2 <- despeckle(hic2, use.scores = 'detrended', focal.size = 2)
hic2
## `HiCExperiment` object with 168,785 contacts over 150 regions
## -------
## fileName: "/home/biocbuild/.cache/R/ExperimentHub/1a8be61e27f712_7752"
## focus: "II:400,000-700,000"
## resolutions(5): 1000 2000 4000 8000 16000
## active resolution: 2000
## interactions: 11325
## scores(5): count balanced expected detrended detrended.despeckled
## topologicalFeatures: compartments(0) borders(0) loops(0) viewpoints(0)
## pairsFile: N/A
## metadata(0):
The added <use.scores>.despeckled
scores correspond to scores averaged using a window, whose width is provided with the focal.size
argument. This results in a smoother Hi-C heatmap, effectively removing the βspecklesβ observed at longer range.
library(InteractionSet)
loops <- system.file('extdata', 'S288C-loops.bedpe', package = 'HiCExperiment') |>
import() |>
makeGInteractionsFromGRangesPairs()
borders <- system.file('extdata', 'S288C-borders.bed', package = 'HiCExperiment') |>
import()
cowplot::plot_grid(
plotMatrix(hic2, caption = FALSE),
plotMatrix(hic2, use.scores = 'detrended', scale = 'linear', limits = c(-1, 1), caption = FALSE),
plotMatrix(
hic2,
use.scores = 'detrended.despeckled',
scale = 'linear',
limits = c(-1, 1),
caption = FALSE,
loops = loops,
borders = borders
),
nrow = 1
)
despeckled
scores
despeckled
scores are in the same scale than the scores
they were computed from.
5.2 Operations between multiple matrices
5.2.1 Merging maps
Hi-C libraries are often sequenced in multiple rounds, for example when high genome coverage is required. This results in multiple contact matrix files being generated. The merge
function can be used to bind several HiCExperiment
objects into a single one.
The different HiCExperiment
objects do not need to all have identical regions
, as shown in the following example.
hic_sub1 <- subsetByOverlaps(hic, GRanges("II:100001-200000"))
hic_sub2 <- subsetByOverlaps(hic, GRanges("II:300001-400000"))
bound_hic <- merge(hic_sub1, hic_sub2)
plotMatrix(bound_hic)
5.2.2 Computing ratio between two maps
Comparing two Hi-C maps can be useful to infer which genomic loci are differentially interacting between experimental conditions. Comparing two HiCExperiment
objects can be done in R
using the divide
function.
For example, we can divide the eco1 mutant Hi-C data by wild-type Hi-C dataset using the divide
function.
hic_eco1 <- import(
CoolFile(HiContactsData('yeast_eco1', 'mcool')),
focus = 'II',
resolution = 2000
)
## see ?HiContactsData and browseVignettes('HiContactsData') for documentation
## loading from cache
div_contacts <- divide(hic_eco1, by = hic)
div_contacts
## `HiCExperiment` object with 996,154 contacts over 407 regions
## -------
## fileName: N/A
## focus: "II"
## resolutions(1): 2000
## active resolution: 2000
## interactions: 60894
## scores(6): count.x balanced.x count.by balanced.by balanced.fc balanced.l2fc
## topologicalFeatures: ()
## pairsFile: N/A
## metadata(2): hce_list operation
We can visually compare wild-type and eco1 maps side by side (left) and their ratio map (right). This highlights the depletion of short-range and increase of long-range interactions in the eco1 dataset.
cowplot::plot_grid(
plotMatrix(hic_eco1, compare.to = hic, limits = c(-4, -1)),
plotMatrix(
div_contacts,
use.scores = 'balanced.fc',
scale = 'log2',
limits = c(-1, 1),
cmap = bwrColors()
)
)
## [1] "/home/biocbuild/.cache/R/ExperimentHub/2dc7284a9c7d2c_7754 | /home/biocbuild/.cache/R/ExperimentHub/1a8be61e27f712_7752"
Session info
sessioninfo::session_info(include_base = TRUE)
## β Session info ββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββ
## setting value
## version R Under development (unstable) (2024-10-21 r87258)
## os Ubuntu 24.04.1 LTS
## system x86_64, linux-gnu
## ui X11
## language (EN)
## collate C
## ctype en_US.UTF-8
## tz America/New_York
## date 2024-11-11
## pandoc 3.1.3 @ /usr/bin/ (via rmarkdown)
##
## β Packages ββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββ
## package * version date (UTC) lib source
## abind 1.4-8 2024-09-12 [2] CRAN (R 4.5.0)
## AnnotationDbi 1.69.0 2024-11-10 [2] Bioconductor 3.21 (R 4.5.0)
## AnnotationHub * 3.15.0 2024-11-10 [2] Bioconductor 3.21 (R 4.5.0)
## backports 1.5.0 2024-05-23 [2] CRAN (R 4.5.0)
## base * 4.5.0 2024-10-23 [3] local
## base64enc 0.1-3 2015-07-28 [2] CRAN (R 4.5.0)
## beeswarm 0.4.0 2021-06-01 [2] CRAN (R 4.5.0)
## Biobase * 2.67.0 2024-11-10 [2] Bioconductor 3.21 (R 4.5.0)
## BiocFileCache * 2.15.0 2024-11-10 [2] Bioconductor 3.21 (R 4.5.0)
## BiocGenerics * 0.53.1 2024-11-10 [2] Bioconductor 3.21 (R 4.5.0)
## BiocIO 1.17.0 2024-11-10 [2] Bioconductor 3.21 (R 4.5.0)
## BiocManager 1.30.25 2024-08-28 [2] CRAN (R 4.5.0)
## BiocParallel 1.41.0 2024-11-10 [2] Bioconductor 3.21 (R 4.5.0)
## BiocVersion 3.21.1 2024-11-10 [2] Bioconductor 3.21 (R 4.5.0)
## Biostrings 2.75.1 2024-11-10 [2] Bioconductor 3.21 (R 4.5.0)
## bit 4.5.0 2024-09-20 [2] CRAN (R 4.5.0)
## bit64 4.5.2 2024-09-22 [2] CRAN (R 4.5.0)
## bitops 1.0-9 2024-10-03 [2] CRAN (R 4.5.0)
## blob 1.2.4 2023-03-17 [2] CRAN (R 4.5.0)
## cachem 1.1.0 2024-05-16 [2] CRAN (R 4.5.0)
## Cairo 1.6-2 2023-11-28 [2] CRAN (R 4.5.0)
## checkmate 2.3.2 2024-07-29 [2] CRAN (R 4.5.0)
## cli 3.6.3 2024-06-21 [2] CRAN (R 4.5.0)
## cluster 2.1.6 2023-12-01 [3] CRAN (R 4.5.0)
## codetools 0.2-20 2024-03-31 [3] CRAN (R 4.5.0)
## colorspace 2.1-1 2024-07-26 [2] CRAN (R 4.5.0)
## compiler 4.5.0 2024-10-23 [3] local
## cowplot 1.1.3 2024-01-22 [2] CRAN (R 4.5.0)
## crayon 1.5.3 2024-06-20 [2] CRAN (R 4.5.0)
## curl 6.0.0 2024-11-05 [2] CRAN (R 4.5.0)
## data.table 1.16.2 2024-10-10 [2] CRAN (R 4.5.0)
## datasets * 4.5.0 2024-10-23 [3] local
## DBI 1.2.3 2024-06-02 [2] CRAN (R 4.5.0)
## dbplyr * 2.5.0 2024-03-19 [2] CRAN (R 4.5.0)
## DelayedArray 0.33.1 2024-11-10 [2] Bioconductor 3.21 (R 4.5.0)
## digest 0.6.37 2024-08-19 [2] CRAN (R 4.5.0)
## doParallel 1.0.17 2022-02-07 [2] CRAN (R 4.5.0)
## dplyr 1.1.4 2023-11-17 [2] CRAN (R 4.5.0)
## dynamicTreeCut 1.63-1 2016-03-11 [2] CRAN (R 4.5.0)
## evaluate 1.0.1 2024-10-10 [2] CRAN (R 4.5.0)
## ExperimentHub * 2.15.0 2024-11-10 [2] Bioconductor 3.21 (R 4.5.0)
## fansi 1.0.6 2023-12-08 [2] CRAN (R 4.5.0)
## farver 2.1.2 2024-05-13 [2] CRAN (R 4.5.0)
## fastcluster 1.2.6 2024-01-12 [2] CRAN (R 4.5.0)
## fastmap 1.2.0 2024-05-15 [2] CRAN (R 4.5.0)
## filelock 1.0.3 2023-12-11 [2] CRAN (R 4.5.0)
## foreach 1.5.2 2022-02-02 [2] CRAN (R 4.5.0)
## foreign 0.8-87 2024-06-26 [3] CRAN (R 4.5.0)
## Formula 1.2-5 2023-02-24 [2] CRAN (R 4.5.0)
## generics * 0.1.3 2022-07-05 [2] CRAN (R 4.5.0)
## GenomeInfoDb * 1.43.0 2024-11-10 [2] Bioconductor 3.21 (R 4.5.0)
## GenomeInfoDbData 1.2.13 2024-10-23 [2] Bioconductor
## GenomicAlignments 1.43.0 2024-11-10 [2] Bioconductor 3.21 (R 4.5.0)
## GenomicRanges * 1.59.0 2024-11-10 [2] Bioconductor 3.21 (R 4.5.0)
## ggbeeswarm 0.7.2 2023-04-29 [2] CRAN (R 4.5.0)
## ggplot2 * 3.5.1 2024-04-23 [2] CRAN (R 4.5.0)
## ggrastr 1.0.2 2023-06-01 [2] CRAN (R 4.5.0)
## glue 1.8.0 2024-09-30 [2] CRAN (R 4.5.0)
## GO.db 3.20.0 2024-10-24 [2] Bioconductor
## graphics * 4.5.0 2024-10-23 [3] local
## grDevices * 4.5.0 2024-10-23 [3] local
## grid 4.5.0 2024-10-23 [3] local
## gridExtra 2.3 2017-09-09 [2] CRAN (R 4.5.0)
## gtable 0.3.6 2024-10-25 [2] CRAN (R 4.5.0)
## HiCExperiment * 1.7.0 2024-11-10 [2] Bioconductor 3.21 (R 4.5.0)
## HiContacts * 1.9.0 2024-11-10 [2] Bioconductor 3.21 (R 4.5.0)
## HiContactsData * 1.9.0 2024-11-07 [2] Bioconductor 3.21 (R 4.5.0)
## Hmisc 5.2-0 2024-10-28 [2] CRAN (R 4.5.0)
## hms 1.1.3 2023-03-21 [2] CRAN (R 4.5.0)
## htmlTable 2.4.3 2024-07-21 [2] CRAN (R 4.5.0)
## htmltools 0.5.8.1 2024-04-04 [2] CRAN (R 4.5.0)
## htmlwidgets 1.6.4 2023-12-06 [2] CRAN (R 4.5.0)
## httr 1.4.7 2023-08-15 [2] CRAN (R 4.5.0)
## impute 1.81.0 2024-11-10 [2] Bioconductor 3.21 (R 4.5.0)
## InteractionSet * 1.35.0 2024-11-10 [2] Bioconductor 3.21 (R 4.5.0)
## IRanges * 2.41.0 2024-11-10 [2] Bioconductor 3.21 (R 4.5.0)
## iterators 1.0.14 2022-02-05 [2] CRAN (R 4.5.0)
## jsonlite 1.8.9 2024-09-20 [2] CRAN (R 4.5.0)
## KEGGREST 1.47.0 2024-11-10 [2] Bioconductor 3.21 (R 4.5.0)
## knitr 1.49 2024-11-08 [2] CRAN (R 4.5.0)
## labeling 0.4.3 2023-08-29 [2] CRAN (R 4.5.0)
## lattice 0.22-6 2024-03-20 [3] CRAN (R 4.5.0)
## lifecycle 1.0.4 2023-11-07 [2] CRAN (R 4.5.0)
## magrittr 2.0.3 2022-03-30 [2] CRAN (R 4.5.0)
## Matrix 1.7-1 2024-10-18 [3] CRAN (R 4.5.0)
## MatrixGenerics * 1.19.0 2024-11-10 [2] Bioconductor 3.21 (R 4.5.0)
## matrixStats * 1.4.1 2024-09-08 [2] CRAN (R 4.5.0)
## memoise 2.0.1 2021-11-26 [2] CRAN (R 4.5.0)
## methods * 4.5.0 2024-10-23 [3] local
## mime 0.12 2021-09-28 [2] CRAN (R 4.5.0)
## munsell 0.5.1 2024-04-01 [2] CRAN (R 4.5.0)
## nnet 7.3-19 2023-05-03 [3] CRAN (R 4.5.0)
## parallel 4.5.0 2024-10-23 [3] local
## pillar 1.9.0 2023-03-22 [2] CRAN (R 4.5.0)
## pkgconfig 2.0.3 2019-09-22 [2] CRAN (R 4.5.0)
## png 0.1-8 2022-11-29 [2] CRAN (R 4.5.0)
## preprocessCore 1.69.0 2024-11-10 [2] Bioconductor 3.21 (R 4.5.0)
## purrr 1.0.2 2023-08-10 [2] CRAN (R 4.5.0)
## R6 2.5.1 2021-08-19 [2] CRAN (R 4.5.0)
## rappdirs 0.3.3 2021-01-31 [2] CRAN (R 4.5.0)
## Rcpp 1.0.13-1 2024-11-02 [2] CRAN (R 4.5.0)
## RCurl 1.98-1.16 2024-07-11 [2] CRAN (R 4.5.0)
## readr 2.1.5 2024-01-10 [2] CRAN (R 4.5.0)
## restfulr 0.0.15 2022-06-16 [2] CRAN (R 4.5.0)
## rhdf5 2.51.0 2024-11-10 [2] Bioconductor 3.21 (R 4.5.0)
## rhdf5filters 1.19.0 2024-11-10 [2] Bioconductor 3.21 (R 4.5.0)
## Rhdf5lib 1.29.0 2024-11-10 [2] Bioconductor 3.21 (R 4.5.0)
## rjson 0.2.23 2024-09-16 [2] CRAN (R 4.5.0)
## rlang 1.1.4 2024-06-04 [2] CRAN (R 4.5.0)
## rmarkdown 2.29 2024-11-04 [2] CRAN (R 4.5.0)
## rpart 4.1.23 2023-12-05 [3] CRAN (R 4.5.0)
## Rsamtools 2.23.0 2024-11-10 [2] Bioconductor 3.21 (R 4.5.0)
## RSpectra 0.16-2 2024-07-18 [2] CRAN (R 4.5.0)
## RSQLite 2.3.7 2024-05-27 [2] CRAN (R 4.5.0)
## rstudioapi 0.17.1 2024-10-22 [2] CRAN (R 4.5.0)
## rtracklayer * 1.67.0 2024-11-10 [2] Bioconductor 3.21 (R 4.5.0)
## S4Arrays 1.7.1 2024-11-10 [2] Bioconductor 3.21 (R 4.5.0)
## S4Vectors * 0.45.1 2024-11-10 [2] Bioconductor 3.21 (R 4.5.0)
## scales 1.3.0 2023-11-28 [2] CRAN (R 4.5.0)
## sessioninfo 1.2.2 2021-12-06 [2] CRAN (R 4.5.0)
## SparseArray 1.7.1 2024-11-10 [2] Bioconductor 3.21 (R 4.5.0)
## splines 4.5.0 2024-10-23 [3] local
## stats * 4.5.0 2024-10-23 [3] local
## stats4 * 4.5.0 2024-10-23 [3] local
## strawr 0.0.92 2024-07-16 [2] CRAN (R 4.5.0)
## stringi 1.8.4 2024-05-06 [2] CRAN (R 4.5.0)
## stringr 1.5.1 2023-11-14 [2] CRAN (R 4.5.0)
## SummarizedExperiment * 1.37.0 2024-11-10 [2] Bioconductor 3.21 (R 4.5.0)
## survival 3.7-0 2024-06-05 [3] CRAN (R 4.5.0)
## terra 1.7-83 2024-10-14 [2] CRAN (R 4.5.0)
## tibble 3.2.1 2023-03-20 [2] CRAN (R 4.5.0)
## tidyr 1.3.1 2024-01-24 [2] CRAN (R 4.5.0)
## tidyselect 1.2.1 2024-03-11 [2] CRAN (R 4.5.0)
## tools 4.5.0 2024-10-23 [3] local
## tzdb 0.4.0 2023-05-12 [2] CRAN (R 4.5.0)
## UCSC.utils 1.3.0 2024-11-10 [2] Bioconductor 3.21 (R 4.5.0)
## utf8 1.2.4 2023-10-22 [2] CRAN (R 4.5.0)
## utils * 4.5.0 2024-10-23 [3] local
## vctrs 0.6.5 2023-12-01 [2] CRAN (R 4.5.0)
## vipor 0.4.7 2023-12-18 [2] CRAN (R 4.5.0)
## vroom 1.6.5 2023-12-05 [2] CRAN (R 4.5.0)
## WGCNA 1.73 2024-09-18 [2] CRAN (R 4.5.0)
## withr 3.0.2 2024-10-28 [2] CRAN (R 4.5.0)
## xfun 0.49 2024-10-31 [2] CRAN (R 4.5.0)
## XML 3.99-0.17 2024-06-25 [2] CRAN (R 4.5.0)
## XVector 0.47.0 2024-11-10 [2] Bioconductor 3.21 (R 4.5.0)
## yaml 2.3.10 2024-07-26 [2] CRAN (R 4.5.0)
## zlibbioc 1.53.0 2024-11-10 [2] Bioconductor 3.21 (R 4.5.0)
##
## [1] /tmp/Rtmp04Rwt0/Rinst31e971dea8c2a
## [2] /home/biocbuild/bbs-3.21-bioc/R/site-library
## [3] /home/biocbuild/bbs-3.21-bioc/R/library
##
## βββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββββ