linkSet 0.99.19
8 August 2025
linkSet 0.99.19
Understanding the 3D architecture of the genome is crucial for deciphering gene regulation and the role of non-coding regions. With techniques like HiC, PCHiC, scATAC-seq, and other high-throughput omics approaches, we can now identify interactions between genes and regulatory elements such as enhancers. These interactions often occur over large genomic distances and can significantly influence gene expression.
The linkSet package addresses a critical need in genomic data analysis: a standardized framework for representing, manipulating, and analyzing gene-enhancer interactions. While several Bioconductor packages deal with genomic interactions, linkSet specifically focuses on gene-enhancer links with specialized methods for this type of data.
Several Bioconductor packages handle genomic interaction data, but linkSet
offers distinct advantages:
InteractionSet: Provides a general representation of genomic interactions but lacks specific tools for gene-enhancer links.
GenomicInteractions: Similar to InteractionSet with visualization capabilities but doesn’t incorporate domain-specific functions for regulatory element analysis.
plyranges: Offers a dplyr-like interface for GRanges objects but isn’t optimized for pairwise genomic interactions.
HiCExperiment: Specializes in Hi-C data but is designed for broader chromatin contacts rather than specific enhancer-promoter interactions.
The linkSet
package fills this gap by providing:
- Specialized methods for annotating promoters
- Functions for calculating interaction metrics relevant to gene regulation
- Distance-based analyses specific to enhancer-promoter interactions
- Visualization tools tailored for regulatory genomic links
linkSet
can seamlessly apply to the following situations:
The linkSet
class is a specialized data structure designed to represent and analyze genomic interactions, particularly focusing on gene-enhancer relationships. It’s built upon proven Bioconductor infrastructure while adding specialized functionality for regulatory genomics.
Key features of the linkSet
class:
The linkSet
class can be constructed from various data types, including GRanges objects for both anchors, or character vectors for bait and GRanges for other ends.
To construct a linkSet
object from GRanges objects:
library(GenomicRanges)
gr1 <- GRanges(
seqnames = c("chr1", "chr1", "chr2"),
ranges = IRanges(start = c(1, 100, 200), width = 10),
strand = "+", symbol = c("Gene1", "Gene2", "Gene3")
)
gr2 <- GRanges(
seqnames = c("chr1", "chr2", "chr2"),
ranges = IRanges(start = c(50, 150, 250), width = 10),
strand = "+"
)
## construct linkSet
linkObj <- linkSet(gr1, gr2, specificCol = "symbol")
linkObj
#> linkSet object with 3 interactions and 1 metadata column:
#> bait seqnames_oe ranges_oe | anchor1.symbol
#> <character> <Rle> <IRanges> | <character>
#> [1] Gene1 --- chr1 50-59 | Gene1
#> [2] Gene2 --- chr2 150-159 | Gene2
#> [3] Gene3 --- chr2 250-259 | Gene3
#> -------
#> regions: 6 ranges and 0 metadata columns
#> seqinfo: 2 sequences from an unspecified genome; no seqlengths
You can also create a linkSet
object from a GInteractions
object using the Convert
function:
library(InteractionSet)
gi <- GInteractions(
anchor1 = c(1, 3, 5),
anchor2 = c(2, 4, 6),
regions = GRanges(
seqnames = c("chr1", "chr1", "chr2", "chr2", "chr3", "chr3"),
ranges = IRanges(start = c(1, 50, 100, 150, 200, 250), width = 10)
)
)
mcols(gi)$score <- c(0.1, 0.2, 0.3)
mcols(gi)$gene <- c("Sox2", "Sox9", "Sox10")
## Convert to linkSet
linkObj_from_gi <- Convert(gi, baitCol = "gene")
linkObj_from_gi
#> linkSet object with 3 interactions and 2 metadata columns:
#> bait seqnames_oe ranges_oe | score gene
#> <character> <Rle> <IRanges> | <numeric> <character>
#> [1] Sox2 --- chr1 50-59 | 0.1 Sox2
#> [2] Sox9 --- chr2 150-159 | 0.2 Sox9
#> [3] Sox10 --- chr3 250-259 | 0.3 Sox10
#> -------
#> regions: 6 ranges and 0 metadata columns
#> seqinfo: 3 sequences from an unspecified genome; no seqlengths
Sometimes, you need to use gene and enhancer information to construct a linkSet
object from GInteractions.
geneGr <- GRanges(
seqnames = c("chr1", "chr2", "chr3"),
ranges = IRanges(start = c(5, 105, 205), width = 20),
symbol = c("Gene1", "Gene2", "Gene3")
)
peakGr <- GRanges(
seqnames = c("chr1", "chr2", "chr3"),
ranges = IRanges(start = c(45, 145, 245), width = 10)
)
# Run baitGInteractions
linkObj_from_gi_2 <- baitGInteractions(gi, geneGr, peakGr, geneSymbol = "symbol")
linkObj_from_gi_2
#> linkSet object with 3 interactions and 3 metadata columns:
#> bait seqnames_oe ranges_oe | anchor1.symbol score
#> <character> <Rle> <IRanges> | <character> <numeric>
#> [1] Gene1 --- chr1 45-54 | Gene1 0.1
#> [2] Gene2 --- chr2 145-154 | Gene2 0.2
#> [3] Gene3 --- chr3 245-254 | Gene3 0.3
#> gene
#> <character>
#> [1] Sox2
#> [2] Sox9
#> [3] Sox10
#> -------
#> regions: 6 ranges and 0 metadata columns
#> seqinfo: 3 sequences from an unspecified genome; no seqlengths
You can also construct a linkSet
object from a data frame.
test_df <- data.frame(
gene = c("gene1", "gene2", "gene3"),
peak = c("chr1:1000-2000", "chr2:3000-4000", "chr3:5000-6000"),
score = c(0.5, 0.7, 0.9),
stringsAsFactors = FALSE
)
# Test successful conversion
linkObj_from_df <- Convert(test_df)
linkObj_from_df
#> linkSet object with 3 interactions and 1 metadata column:
#> bait seqnames_oe ranges_oe | score
#> <character> <Rle> <IRanges> | <numeric>
#> [1] gene1 --- chr1 1000-2000 | 0.5
#> [2] gene2 --- chr2 3000-4000 | 0.7
#> [3] gene3 --- chr3 5000-6000 | 0.9
#> -------
#> regions: 3 ranges and 0 metadata columns
#> seqinfo: 3 sequences from an unspecified genome; no seqlengths
The important elements of the linkSet
object can be easily accessed by the following accessors:
regions(x)
: Get the regions of the linkSet
object.bait(x)
: Get the bait of the linkSet
object.oe(x)
: Get the other end of the linkSet
object.regionsBait(x)
: Get the bait regions of the linkSet
object.regionsOe(x)
: Get the other end regions of the linkSet
object.mcols(x)
: Get the metadata of the linkSet
object.linkSet::regions(linkObj)
#> GRanges object with 6 ranges and 0 metadata columns:
#> seqnames ranges strand
#> <Rle> <IRanges> <Rle>
#> [1] chr1 1-10 +
#> [2] chr1 50-59 +
#> [3] chr1 100-109 +
#> [4] chr2 150-159 +
#> [5] chr2 200-209 +
#> [6] chr2 250-259 +
#> -------
#> seqinfo: 2 sequences from an unspecified genome; no seqlengths
regionsBait(linkObj)
#> GRanges object with 3 ranges and 0 metadata columns:
#> seqnames ranges strand
#> <Rle> <IRanges> <Rle>
#> [1] chr1 1-10 +
#> [2] chr1 100-109 +
#> [3] chr2 200-209 +
#> -------
#> seqinfo: 2 sequences from an unspecified genome; no seqlengths
Wait… You set wrong bait or oe region? You can also change it easily:
new_bait <- c("Gene40", "Gene41", "Gene42")
new_oe <- GRanges(
seqnames = c("chr1", "chr2", "chr3"),
ranges = IRanges(start = c(5, 105, 205), width = 20),
symbol = c("Gene1", "Gene2", "Gene3")
)
bait(linkObj) <- new_bait
oe(linkObj) <- new_oe
linkObj
#> linkSet object with 3 interactions and 2 metadata columns:
#> bait seqnames_oe ranges_oe | anchor1.symbol symbol
#> <character> <Rle> <IRanges> | <character> <character>
#> [1] Gene40 --- chr1 5-24 | Gene1 Gene1
#> [2] Gene41 --- chr2 105-124 | Gene2 Gene2
#> [3] Gene42 --- chr3 205-224 | Gene3 Gene3
#> -------
#> regions: 6 ranges and 0 metadata columns
#> seqinfo: 3 sequences from an unspecified genome; no seqlengths
LinkSet object can be easily subsetted by index.
linkset_sub <- linkObj[1:2]
linkset_sub
#> linkSet object with 2 interactions and 2 metadata columns:
#> bait seqnames_oe ranges_oe | anchor1.symbol symbol
#> <character> <Rle> <IRanges> | <character> <character>
#> [1] Gene40 --- chr1 5-24 | Gene1 Gene1
#> [2] Gene41 --- chr2 105-124 | Gene2 Gene2
#> -------
#> regions: 6 ranges and 0 metadata columns
#> seqinfo: 3 sequences from an unspecified genome; no seqlengths
If you are interested in a specific gene or region, you can subset the linkSet
object by the bait or oe region.
linkset_sub_2 <- subsetBait(linkObj, "Gene1")
linkset_sub_2
#> linkSet object with 0 interactions and 2 metadata columns:
#> bait seqnames_oe ranges_oe | anchor1.symbol symbol
#> <character> <Rle> <IRanges> | <character> <character>
#> -------
#> regions: 0 ranges and 0 metadata columns
#> seqinfo: 3 sequences from an unspecified genome; no seqlengths
linkset_sub_3 <- subsetBaitRegion(linkObj, "chr1:100-200")
linkset_sub_3
#> linkSet object with 1 interaction and 2 metadata columns:
#> bait seqnames_oe ranges_oe | anchor1.symbol symbol
#> <character> <Rle> <IRanges> | <character> <character>
#> [1] Gene41 --- chr2 105-124 | Gene2 Gene2
#> -------
#> regions: 2 ranges and 0 metadata columns
#> seqinfo: 3 sequences from an unspecified genome; no seqlengths
You can also concatenate two linkSet
objects.
linkset_concat <- c(linkObj, linkObj)
linkset_concat
#> linkSet object with 6 interactions and 2 metadata columns:
#> bait seqnames_oe ranges_oe | anchor1.symbol symbol
#> <character> <Rle> <IRanges> | <character> <character>
#> [1] Gene40 --- chr1 5-24 | Gene1 Gene1
#> [2] Gene41 --- chr2 105-124 | Gene2 Gene2
#> [3] Gene42 --- chr3 205-224 | Gene3 Gene3
#> [4] Gene40 --- chr1 5-24 | Gene1 Gene1
#> [5] Gene41 --- chr2 105-124 | Gene2 Gene2
#> [6] Gene42 --- chr3 205-224 | Gene3 Gene3
#> -------
#> regions: 6 ranges and 0 metadata columns
#> seqinfo: 3 sequences from an unspecified genome; no seqlengths
The features of linkSet is inherited from GRanges, which means you can use all the functions in GRanges to manipulate the linkSet object. The most wonderful things is that you can change the bait region and oe region separately.
data(linkExample)
linkExample
#> linkSet object with 5 interactions and 1 metadata column:
#> bait seqnames_oe ranges_oe | anchor1.symbol
#> <character> <Rle> <IRanges> | <character>
#> [1] Gene1 --- chr1 50-59 | Gene1
#> [2] Gene1 --- chr2 150-159 | Gene1
#> [3] Gene2 --- chr2 250-259 | Gene2
#> [4] Gene3 --- chr4 350-359 | Gene3
#> [5] Gene3 --- chr4 450-459 | Gene3
#> -------
#> regions: 10 ranges and 0 metadata columns
#> seqinfo: 4 sequences from an unspecified genome; no seqlengths
## resize the bait region
resize_bait <- resizeRegions(linkExample, width = 75, fix = "start", region = "bait")
resize_bait
#> linkSet object with 5 interactions and 1 metadata column:
#> bait seqnames_oe ranges_oe | anchor1.symbol
#> <character> <Rle> <IRanges> | <character>
#> [1] Gene1 --- chr1 50-59 | Gene1
#> [2] Gene1 --- chr2 150-159 | Gene1
#> [3] Gene2 --- chr2 250-259 | Gene2
#> [4] Gene3 --- chr4 350-359 | Gene3
#> [5] Gene3 --- chr4 450-459 | Gene3
#> -------
#> regions: 10 ranges and 0 metadata columns
#> seqinfo: 4 sequences from an unspecified genome; no seqlengths
## resize the oe region
resize_oe <- resizeRegions(linkExample, width = 75, fix = "start", region = "oe")
resize_oe
#> linkSet object with 5 interactions and 1 metadata column:
#> bait seqnames_oe ranges_oe | anchor1.symbol
#> <character> <Rle> <IRanges> | <character>
#> [1] Gene1 --- chr1 50-124 | Gene1
#> [2] Gene1 --- chr2 150-224 | Gene1
#> [3] Gene2 --- chr2 250-324 | Gene2
#> [4] Gene3 --- chr4 350-424 | Gene3
#> [5] Gene3 --- chr4 450-524 | Gene3
#> -------
#> regions: 10 ranges and 0 metadata columns
#> seqinfo: 4 sequences from an unspecified genome; no seqlengths
The reduce
method allows the linkSet
objects to be collapsed across the whole of the linkSet
object.
reduce_linkset <- reduce(linkset_concat)
reduce_linkset
#> linkSet object with 3 interactions and 3 metadata columns:
#> bait seqnames_oe ranges_oe | anchor1.symbol symbol
#> <character> <Rle> <IRanges> | <character> <character>
#> [1] Gene40 --- chr1 5-24 | Gene1 Gene1
#> [2] Gene41 --- chr2 105-124 | Gene2 Gene2
#> [3] Gene42 --- chr3 205-224 | Gene3 Gene3
#> count
#> <integer>
#> [1] 2
#> [2] 2
#> [3] 2
#> -------
#> regions: 6 ranges and 0 metadata columns
#> seqinfo: 3 sequences from an unspecified genome; no seqlengths
By default, we will collapse the duplicate interactions, and count the interactions and store in the metadata columns.
You can also choose not to count the interactions by set countInteractions = FALSE
.
reduce_linkset_2 <- linkSet::reduceRegions(linkset_concat, countInteractions = FALSE)
reduce_linkset_2
#> linkSet object with 6 interactions and 2 metadata columns:
#> bait seqnames_oe ranges_oe | anchor1.symbol symbol
#> <character> <Rle> <IRanges> | <character> <character>
#> [1] Gene40 --- chr1 5-24 | Gene1 Gene1
#> [2] Gene41 --- chr2 105-124 | Gene2 Gene2
#> [3] Gene42 --- chr3 205-224 | Gene3 Gene3
#> [4] Gene40 --- chr1 5-24 | Gene1 Gene1
#> [5] Gene41 --- chr2 105-124 | Gene2 Gene2
#> [6] Gene42 --- chr3 205-224 | Gene3 Gene3
#> -------
#> regions: 6 ranges and 0 metadata columns
#> seqinfo: 3 sequences from an unspecified genome; no seqlengths
There are two metrics to diagnose the quality of the linkSet
object:
1. The intra/inter-chromosomal interactions.
2. The distance between the bait and oe region.
The diagnose
function can help you to diagnose the linkSet
object.
#> linkSet object with 5 interactions and 3 metadata columns:
#> bait seqnames_oe ranges_oe | anchor1.symbol inter_type
#> <character> <Rle> <IRanges> | <character> <character>
#> [1] Gene1 --- chr1 50-59 | Gene1 inter
#> [2] Gene1 --- chr2 150-159 | Gene1 intra
#> [3] Gene2 --- chr2 250-259 | Gene2 inter
#> [4] Gene3 --- chr4 350-359 | Gene3 intra
#> [5] Gene3 --- chr4 450-459 | Gene3 inter
#> distance
#> <integer>
#> [1] 49
#> [2] <NA>
#> [3] 50
#> [4] <NA>
#> [5] 50
#> -------
#> regions: 10 ranges and 0 metadata columns
#> seqinfo: 4 sequences from an unspecified genome; no seqlengths
linkset_concat
#> linkSet object with 6 interactions and 2 metadata columns:
#> bait seqnames_oe ranges_oe | anchor1.symbol symbol
#> <character> <Rle> <IRanges> | <character> <character>
#> [1] Gene40 --- chr1 5-24 | Gene1 Gene1
#> [2] Gene41 --- chr2 105-124 | Gene2 Gene2
#> [3] Gene42 --- chr3 205-224 | Gene3 Gene3
#> [4] Gene40 --- chr1 5-24 | Gene1 Gene1
#> [5] Gene41 --- chr2 105-124 | Gene2 Gene2
#> [6] Gene42 --- chr3 205-224 | Gene3 Gene3
#> -------
#> regions: 6 ranges and 0 metadata columns
#> seqinfo: 3 sequences from an unspecified genome; no seqlengths
The annotated distance and interaction type is shown in the distance
and inter_type
column.
You can remove the intrachromosomal interactions by:
linkset_concat <- filterLinks(linkset_concat, filter_intra = TRUE)
linkset_concat
#> linkSet object with 2 interactions and 4 metadata columns:
#> bait seqnames_oe ranges_oe | anchor1.symbol symbol
#> <character> <Rle> <IRanges> | <character> <character>
#> [1] Gene40 --- chr1 5-24 | Gene1 Gene1
#> [2] Gene40 --- chr1 5-24 | Gene1 Gene1
#> inter_type distance
#> <character> <integer>
#> [1] inter 9
#> [2] inter 9
#> -------
#> regions: 2 ranges and 0 metadata columns
#> seqinfo: 3 sequences from an unspecified genome; no seqlengths
It might be common that you only have the bait name, but don’t have the exact bait region. The annotation function can help you to annotate the bait region from a gene annotation.
gr1 <- GRanges(
seqnames = c("chr1", "chr3", "chr3"),
ranges = IRanges(start = c(1000, 2000, 3000), width = 100),
strand = "+", symbol = c("BRCA1", "TP53", "NONEXISTENT")
)
gr2 <- GRanges(
seqnames = c("chr1", "chr2", "chr3"),
ranges = IRanges(start = c(5000, 6000, 7000), width = 100),
strand = "+"
)
linkObj <- linkSet(gr1, gr2, specificCol = "symbol")
# Test annotatePromoter
annotated_linkset <- suppressWarnings(annotatePromoter(linkObj, genome = "hg38", upstream = 500, overwrite = TRUE))
annotated_linkset
#> linkSet object with 3 interactions and 1 metadata column:
#> bait seqnames_oe ranges_oe | anchor1.symbol
#> <character> <Rle> <IRanges> | <character>
#> [1] BRCA1 --- chr1 5000-5099 | BRCA1
#> [2] TP53 --- chr2 6000-6099 | TP53
#> [3] NONEXISTENT --- chr3 7000-7099 | NONEXISTENT
#> -------
#> regions: 6 ranges and 0 metadata columns
#> seqinfo: 712 sequences (1 circular) from 2 genomes (hg38, NA)
Functional enhancers usually regulate multiple genes, we can use cross gene analysis to identify the cross gene enhancers.
annotated_linkset <- crossGeneEnhancer(annotated_linkset)
annotated_linkset <- orderLinks(annotated_linkset, by = "crossFreq", decreasing = TRUE)
annotated_linkset
#> linkSet object with 3 interactions and 2 metadata columns:
#> bait seqnames_oe ranges_oe | anchor1.symbol crossFreq
#> <character> <Rle> <IRanges> | <character> <integer>
#> [1] BRCA1 --- chr1 5000-5099 | BRCA1 1
#> [2] TP53 --- chr2 6000-6099 | TP53 1
#> [3] NONEXISTENT --- chr3 7000-7099 | NONEXISTENT 1
#> -------
#> regions: 6 ranges and 0 metadata columns
#> seqinfo: 712 sequences (1 circular) from 2 genomes (hg38, NA)
linkSet also implement the CHiCANE
analysis, which can identifying the high-confidence interactions.
CHiCANE needs to count interactibility and dist before running the statistical test.
annotated_linkset <- countInteractibility(annotated_linkset)
annotated_linkset <- linkSet::pairdist(annotated_linkset)
annotated_linkset <- run_chicane(annotated_linkset)
#> [1] "Not found column 'bait.id', adding bait name as default."
#> [1] "Not found column 'bait.to.bait' Set 'FALSE' as default."
#> [1] "fitting model...."
annotated_linkset
#> linkSet object with 3 interactions and 12 metadata columns:
#> bait seqnames_oe ranges_oe | anchor1.symbol crossFreq
#> <character> <Rle> <IRanges> | <character> <integer>
#> [1] BRCA1 --- chr1 5000-5099 | BRCA1 1
#> [2] TP53 --- chr2 6000-6099 | TP53 1
#> [3] NONEXISTENT --- chr3 7000-7099 | NONEXISTENT 1
#> inter_type count bait.trans.count target.trans.count distance
#> <character> <integer> <numeric> <numeric> <integer>
#> [1] intra 1 1 1 <NA>
#> [2] intra 1 1 1 <NA>
#> [3] intra 1 1 1 <NA>
#> bait.id bait.to.bait expected p.value q.value
#> <character> <logical> <numeric> <numeric> <numeric>
#> [1] BRCA1 FALSE 1 0.632121 0.632121
#> [2] TP53 FALSE 1 0.632121 0.632121
#> [3] NONEXISTENT FALSE 1 0.632121 0.632121
#> -------
#> regions: 6 ranges and 0 metadata columns
#> seqinfo: 712 sequences (1 circular) from 2 genomes (hg38, NA)
With plot_genomic_ranges
function, you can visualize the link between bait and oe region.
You can choose the oe-centric views to visualize the functional enhancer which regulate multiple genes.
# Note: visualization requires annotated bait regions
# linkExample doesn't have regionsBait annotated, so this will be skipped
plot_genomic_ranges(linkExample, extend.base = 10)
You can adjust the range of the views.
sessionInfo()
#> R version 4.5.1 (2025-06-13)
#> Platform: x86_64-pc-linux-gnu
#> Running under: Ubuntu 24.04.2 LTS
#>
#> Matrix products: default
#> BLAS: /home/biocbuild/bbs-3.22-bioc/R/lib/libRblas.so
#> LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.12.0 LAPACK version 3.12.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] data.table_1.17.8
#> [2] MASS_7.3-65
#> [3] rtracklayer_1.69.1
#> [4] InteractionSet_1.37.1
#> [5] SummarizedExperiment_1.39.1
#> [6] MatrixGenerics_1.21.0
#> [7] matrixStats_1.5.0
#> [8] Organism.dplyr_1.37.1
#> [9] AnnotationFilter_1.33.0
#> [10] dplyr_1.1.4
#> [11] org.Mm.eg.db_3.21.0
#> [12] TxDb.Mmusculus.UCSC.mm10.knownGene_3.10.0
#> [13] GenomicFeatures_1.61.6
#> [14] AnnotationDbi_1.71.1
#> [15] Biobase_2.69.0
#> [16] linkSet_0.99.19
#> [17] GenomicRanges_1.61.1
#> [18] Seqinfo_0.99.2
#> [19] IRanges_2.43.0
#> [20] S4Vectors_0.47.0
#> [21] BiocGenerics_0.55.1
#> [22] generics_0.1.4
#> [23] BiocStyle_2.37.1
#>
#> loaded via a namespace (and not attached):
#> [1] DBI_1.2.3
#> [2] bitops_1.0-9
#> [3] httr2_1.2.1
#> [4] rlang_1.1.6
#> [5] magrittr_2.0.3
#> [6] compiler_4.5.1
#> [7] RSQLite_2.4.2
#> [8] png_0.1-8
#> [9] vctrs_0.6.5
#> [10] pkgconfig_2.0.3
#> [11] crayon_1.5.3
#> [12] fastmap_1.2.0
#> [13] dbplyr_2.5.0
#> [14] magick_2.8.7
#> [15] XVector_0.49.0
#> [16] labeling_0.4.3
#> [17] Rsamtools_2.25.2
#> [18] rmarkdown_2.29
#> [19] UCSC.utils_1.5.0
#> [20] tinytex_0.57
#> [21] purrr_1.1.0
#> [22] bit_4.6.0
#> [23] xfun_0.52
#> [24] cachem_1.1.0
#> [25] GenomeInfoDb_1.45.9
#> [26] jsonlite_2.0.0
#> [27] blob_1.2.4
#> [28] DelayedArray_0.35.2
#> [29] BiocParallel_1.43.4
#> [30] parallel_4.5.1
#> [31] R6_2.6.1
#> [32] bslib_0.9.0
#> [33] RColorBrewer_1.1-3
#> [34] jquerylib_0.1.4
#> [35] Rcpp_1.1.0
#> [36] bookdown_0.43
#> [37] iterators_1.0.14
#> [38] knitr_1.50
#> [39] Matrix_1.7-3
#> [40] tidyselect_1.2.1
#> [41] dichromat_2.0-0.1
#> [42] abind_1.4-8
#> [43] yaml_2.3.10
#> [44] codetools_0.2-20
#> [45] curl_6.4.0
#> [46] lattice_0.22-7
#> [47] tibble_3.3.0
#> [48] withr_3.0.2
#> [49] KEGGREST_1.49.1
#> [50] evaluate_1.0.4
#> [51] BiocFileCache_2.99.5
#> [52] Biostrings_2.77.2
#> [53] pillar_1.11.0
#> [54] BiocManager_1.30.26
#> [55] filelock_1.0.3
#> [56] foreach_1.5.2
#> [57] RCurl_1.98-1.17
#> [58] ggplot2_3.5.2
#> [59] scales_1.4.0
#> [60] glue_1.8.0
#> [61] lazyeval_0.2.2
#> [62] tools_4.5.1
#> [63] BiocIO_1.19.0
#> [64] TxDb.Hsapiens.UCSC.hg38.knownGene_3.21.0
#> [65] GenomicAlignments_1.45.2
#> [66] XML_3.99-0.18
#> [67] grid_4.5.1
#> [68] patchwork_1.3.1
#> [69] restfulr_0.0.16
#> [70] cli_3.6.5
#> [71] rappdirs_0.3.3
#> [72] S4Arrays_1.9.1
#> [73] gtable_0.3.6
#> [74] sass_0.4.10
#> [75] digest_0.6.37
#> [76] SparseArray_1.9.1
#> [77] org.Hs.eg.db_3.21.0
#> [78] rjson_0.2.23
#> [79] farver_2.1.2
#> [80] memoise_2.0.1
#> [81] htmltools_0.5.8.1
#> [82] lifecycle_1.0.4
#> [83] httr_1.4.7
#> [84] bit64_4.6.0-1