SGSeq provides functions and data structures for analyzing splice events from RNA-seq data. An overview of SGSeq classes and how they are related is shown in Figure 1 and summarized below:
If you use SGSeq, please cite:
When starting a new project, SGSeq requires information about the samples to be analyzed. This information is obtained once initially, and can then be used for all future analyses. Sample information is provided as a data frame with the following columns:
Sample information can be stored in a data.frame or DataFrame object (if BAM files are specified as a BamFileList, it must be stored in a DataFrame). Sample information can be obtained automatically with function getBamInfo(), which takes as input a data frame with columns sample_name and file_bam and extracts all other information from the specified BAM files.
For SGSeq to work correctly it is essential that reads were mapped with a splice-aware alignment program, such as GSNAP (Wu and Nacu 2010), HISAT (Kim, Langmead, and Salzberg 2015) or STAR (Dobin et al. 2013), which generates SAM/BAM files with a custom tag ‘XS’ for spliced reads, indicating the direction of transcription. BAM files must be indexed. Note that lib_size should be the total number of aligned fragments, even if BAM files were subset to genomic regions of interest. The total number of fragments is required for accurate calculation of FPKM values (fragments per kilobase and million sequenced fragments). Here, the term ‘fragment’ denotes a sequenced cDNA fragment, which is represented by a single read in single-end data, or a pair of reads in paired-end data.
This vignette illustrates an analysis of paired-end RNA-seq data from four tumor and four normal colorectal samples, which are part of a data set published in (Seshagiri et al. 2012). RNA-seq reads were mapped to the human reference genome using GSNAP (Wu and Nacu 2010). The analysis is based on BAM files that were subset to reads mapping to a single gene of interest (FBXO31). A data.frame si with sample information was generated from the original BAM files with function getBamInfo(). Note that column lib_size reflects the total number of aligned fragments in the original BAM files.
## sample_name file_bam paired_end read_length frag_length lib_size
## 1 N1 N1.bam TRUE 75 293 12405197
## 2 N2 N2.bam TRUE 75 197 13090179
## 3 N3 N3.bam TRUE 75 206 14983084
## 4 N4 N4.bam TRUE 75 207 15794088
## 5 T1 T1.bam TRUE 75 284 14345976
## 6 T2 T2.bam TRUE 75 235 15464168
## 7 T3 T3.bam TRUE 75 259 15485954
## 8 T4 T4.bam TRUE 75 247 15808356
The following code block sets the correct BAM file paths for the current SGSeq installation.
Transcript annotation can be obtained from a TxDb object or imported from GFF format using function importTranscripts(). Alternatively, transcripts can be specified as a GRangesList of exons grouped by transcripts. In the following code block, the UCSC knownGene table is loaded as a TxDb object. Transcripts on chromosome 16 (where the FBXO31 gene is located) are retained, and chromosome names are changed to match the naming convention used in the BAM files.
library(TxDb.Hsapiens.UCSC.hg19.knownGene)
txdb <- TxDb.Hsapiens.UCSC.hg19.knownGene
txdb <- keepSeqlevels(txdb, "chr16")
seqlevelsStyle(txdb) <- "NCBI"
To work with the annotation in the SGSeq framework, transcript features are extracted from the TxDb object using function convertToTxFeatures(). Only features overlapping the genomic locus of the FBXO31 gene are retained. The genomic coordinates of FBXO31 are stored in a GRanges object gr.
## TxFeatures object with 6 ranges and 0 metadata columns:
## seqnames ranges strand type
## <Rle> <IRanges> <Rle> <factor>
## [1] 16 87362942-87365116 - L
## [2] 16 87365116-87367492 - J
## [3] 16 87367492-87367892 - I
## [4] 16 87367892-87368910 - J
## [5] 16 87368910-87369063 - I
## [6] 16 87369063-87369761 - J
## txName geneName
## <CharacterList> <CharacterList>
## [1] uc002fjv.3,uc002fjw.3,uc010vot.2 79791
## [2] uc002fjv.3,uc002fjw.3,uc010vot.2 79791
## [3] uc002fjv.3,uc002fjw.3,uc010vot.2 79791
## [4] uc002fjv.3,uc002fjw.3,uc010vot.2 79791
## [5] uc002fjv.3,uc002fjw.3,uc010vot.2 79791
## [6] uc002fjv.3,uc002fjw.3,uc010vot.2 79791
## -------
## seqinfo: 1 sequence from GRCh37.p13 genome
convertToTxFeatures() returns a TxFeatures object, which is a GRanges-like object with additional columns. Column type indicates the feature type and can take values
Columns txName and geneName indicate the transcript and gene each feature belongs to. Note that a feature can belong to more than one transcript, and accordingly these columns can store multiple values for each feature. For TxFeatures and other data structures in SGSeq, columns can be accessed with accessor functions as shown in the following code block.
## [1] L J I J I J I J I J I J I J I F J J F U I J F
## Levels: J I F L U
## CharacterList of length 6
## [[1]] uc002fjv.3 uc002fjw.3 uc010vot.2
## [[2]] uc002fjv.3 uc002fjw.3 uc010vot.2
## [[3]] uc002fjv.3 uc002fjw.3 uc010vot.2
## [[4]] uc002fjv.3 uc002fjw.3 uc010vot.2
## [[5]] uc002fjv.3 uc002fjw.3 uc010vot.2
## [[6]] uc002fjv.3 uc002fjw.3 uc010vot.2
## CharacterList of length 6
## [[1]] 79791
## [[2]] 79791
## [[3]] 79791
## [[4]] 79791
## [[5]] 79791
## [[6]] 79791
Exons stored in a TxFeatures object represent the exons spliced together in an RNA molecule. In the context of the splice graph, exons are represented by unique non-overlapping exonic regions. Function convertToSGFeatures() converts TxFeatures to SGFeatures. In the process, overlapping exons are disjoined into non-overlapping exon bins.
## SGFeatures object with 6 ranges and 0 metadata columns:
## seqnames ranges strand type splice5p splice3p featureID
## <Rle> <IRanges> <Rle> <factor> <logical> <logical> <integer>
## [1] 16 87362942-87365116 - E TRUE FALSE 1
## [2] 16 87365116 - A <NA> <NA> 2
## [3] 16 87365116-87367492 - J <NA> <NA> 3
## [4] 16 87367492 - D <NA> <NA> 4
## [5] 16 87367492-87367892 - E TRUE TRUE 5
## [6] 16 87367892 - A <NA> <NA> 6
## geneID txName geneName
## <integer> <CharacterList> <CharacterList>
## [1] 1 uc002fjv.3,uc002fjw.3,uc010vot.2 79791
## [2] 1 uc002fjv.3,uc002fjw.3,uc010vot.2 79791
## [3] 1 uc002fjv.3,uc002fjw.3,uc010vot.2 79791
## [4] 1 uc002fjv.3,uc002fjw.3,uc010vot.2 79791
## [5] 1 uc002fjv.3,uc002fjw.3,uc010vot.2 79791
## [6] 1 uc002fjv.3,uc002fjw.3,uc010vot.2 79791
## -------
## seqinfo: 1 sequence from GRCh37.p13 genome
Similar to TxFeatures, SGFeatures are GRanges-like objects with additional columns. Column type for an SGFeatures object takes values
By convention, splice donor and acceptor sites correspond to the exonic positions immediately flanking the intron. SGFeatures has additional columns not included in TxFeatures: spliced5p and spliced3p indicate whether exon bins have a mandatory splice at the 5′ and 3′ end, respectively. This information is used to determine whether a read is structurally compatible with an exon bin, and whether an exon bin is consistent with an annotated transcript. featureID provides a unique identifier for each feature, geneID indicates the unique component of the splice graph a feature belongs to.
This section illustrates an analysis based on annotated transcripts. Function analyzeFeatures() converts transcript features to splice graph features and obtains compatible fragment counts for each feature and each sample.
## class: SGFeatureCounts
## dim: 42 8
## metadata(0):
## assays(2): counts FPKM
## rownames: NULL
## rowData names(0):
## colnames(8): N1 N2 ... T3 T4
## colData names(6): sample_name file_bam ... frag_length lib_size
analyzeFeatures() returns an SGFeatureCounts object. SGFeatureCounts contains the sample information as colData, splice graph features as rowRanges and assays counts and FPKM, which store compatible fragment counts and FPKMs, respectively. The different data types can be accessed using accessor functions with the same name.
Counts for exons and splice junctions are based on structurally compatible fragments. In the case of splice donors and acceptors, counts indicate the number of fragments with reads spanning the spliced boundary (overlapping the splice site and the flanking intronic position).
FPKM values are calculated as $\frac{x}{NL}10^6$, where x is the number of compatible fragments, N is the library size (stored in lib_size) and L is the effective feature length (the number of possible positions for a compatible fragment). For paired-end data it is assumed that fragment length is equal to frag_length.
FPKMs for splice graph features can be visualized with function plotFeatures. plotFeatures generates a two-panel figure with a splice graph shown in the top panel and a heatmap of expression levels for individual features in the bottom panel. For customization of plotFeatures output, see section Visualization. The plotting function invisibly returns a data.frame with details about the splice graph shown in the plot.
Note that the splice graph includes three alternative transcript start sites (TSSs). However, the heatmap indicates that the first TSS is not used in this data set.
Instead of relying on existing annotation, annotation can be augmented with predictions from RNA-seq data, or the splice graph can be constructed from RNA-seq data without the use of annotation. The following code block predicts transcript features supported by RNA-seq reads, converts them into splice graph features, and then obtains compatible fragment counts. For details on how predictions are obtained, please see (Goldstein et al. 2016).
## SGFeatures object with 6 ranges and 0 metadata columns:
## seqnames ranges strand type splice5p splice3p featureID
## <Rle> <IRanges> <Rle> <factor> <logical> <logical> <integer>
## [1] 16 87362930-87365116 - E TRUE FALSE 1
## [2] 16 87365116 - A <NA> <NA> 2
## [3] 16 87365116-87367492 - J <NA> <NA> 3
## [4] 16 87367492 - D <NA> <NA> 4
## [5] 16 87367492-87367892 - E TRUE TRUE 5
## [6] 16 87367892 - A <NA> <NA> 6
## geneID txName geneName
## <integer> <CharacterList> <CharacterList>
## [1] 1
## [2] 1
## [3] 1
## [4] 1
## [5] 1
## [6] 1
## -------
## seqinfo: 84 sequences from an unspecified genome
For interpretation, predicted features can be annotated with respect to known transcripts. The annotate() function assigns compatible transcripts to each feature and stores the corresponding transcript and gene name in columns txName and geneName, respectively.
## SGFeatures object with 6 ranges and 0 metadata columns:
## seqnames ranges strand type splice5p splice3p featureID
## <Rle> <IRanges> <Rle> <factor> <logical> <logical> <integer>
## [1] 16 87362930-87365116 - E TRUE FALSE 1
## [2] 16 87365116 - A <NA> <NA> 2
## [3] 16 87365116-87367492 - J <NA> <NA> 3
## [4] 16 87367492 - D <NA> <NA> 4
## [5] 16 87367492-87367892 - E TRUE TRUE 5
## [6] 16 87367892 - A <NA> <NA> 6
## geneID txName geneName
## <integer> <CharacterList> <CharacterList>
## [1] 1 uc002fjv.3,uc002fjw.3,uc010vot.2 79791
## [2] 1 uc002fjv.3,uc002fjw.3,uc010vot.2 79791
## [3] 1 uc002fjv.3,uc002fjw.3,uc010vot.2 79791
## [4] 1 uc002fjv.3,uc002fjw.3,uc010vot.2 79791
## [5] 1 uc002fjv.3,uc002fjw.3,uc010vot.2 79791
## [6] 1 uc002fjv.3,uc002fjw.3,uc010vot.2 79791
## -------
## seqinfo: 84 sequences from an unspecified genome
The predicted splice graph and FPKMs can be visualized as previously. Splice graph features with missing annotation are highlighted using argument color_novel.
Note that most exons and splice junctions predicted from RNA-seq data are consistent with transcripts in the UCSC knownGene table (shown in gray). However, in contrast to the previous figure, the predicted gene model does not include parts of the splice graph that are not expressed in the data. Also, an unannotated exon (E3, shown in red) was discovered from the RNA-seq data, which is expressed in three of the four normal colorectal samples (N2, N3, N4).
Instead of considering the full splice graph of a gene, the analysis can be focused on individual splice events. Function analyzeVariants() recursively identifies splice events from the graph, obtains representative counts for each splice variant, and computes estimates of relative splice variant usage, also referred to as ‘percentage spliced in’ (PSI or Ψ) Katz et al. (2010).
## class: SGVariantCounts
## dim: 2 8
## metadata(0):
## assays(5): countsVariant5p countsVariant3p countsEvent5p countsEvent3p
## variantFreq
## rownames: NULL
## rowData names(20): from to ... variantType variantName
## colnames(8): N1 N2 ... T3 T4
## colData names(6): sample_name file_bam ... frag_length lib_size
analyzeVariants() returns an SGVariantCounts object. Sample information is stored as colData, and SGVariants as rowRanges. Assay variantFreq stores estimates of relative usage for each splice variant and sample. As previously, the different data types can be accessed using accessor functions. Information on splice variants is stored in SGVariants metadata columns and can be accessed with function mcols() as shown below. For a detailed description of columns, see the manual page for SGVariants.
## DataFrame with 2 rows and 20 columns
## from to type featureID segmentID closed5p
## <character> <character> <character> <character> <character> <logical>
## 1 D:16:87393901:- A:16:87380856:- J 28 4 TRUE
## 2 D:16:87393901:- A:16:87380856:- JEJ 32,30,27 2 TRUE
## closed3p closed5pEvent closed3pEvent geneID eventID variantID
## <logical> <logical> <logical> <integer> <integer> <integer>
## 1 TRUE TRUE TRUE 1 1 1
## 2 TRUE TRUE TRUE 1 1 2
## featureID5p featureID3p featureID5pEvent featureID3pEvent
## <IntegerList> <IntegerList> <IntegerList> <IntegerList>
## 1 28 28 28,32 28,27
## 2 32 27 28,32 28,27
## txName geneName variantType
## <CharacterList> <CharacterList> <CharacterList>
## 1 uc002fjv.3,uc002fjw.3,uc010vot.2 79791 SE:S
## 2 79791 SE:I
## variantName
## <character>
## 1 79791_1_1/2_SE
## 2 79791_1_2/2_SE
Splice variants are quantified locally, based on structurally compatible fragments that extend across the start or end of each variant. Local estimates of relative usage Ψi for variant i are obtained as the number of fragments compatible with i divided by the number of fragments compatible with any variant belonging to the same event. For variant start S and variant end E, Ψ̂iS = xiS/x.S and Ψ̂iE = xiE/x.E, respectively. For variants with valid estimates Ψ̂iS and Ψ̂iE, a single estimate is calculated as a weighted mean of local estimates Ψ̂i = x.S/(x.S + x.E)Ψ̂iS + x.E/(x.S + x.E)Ψ̂iE.
Estimates of relative usage can be unreliable for events with low read count. If argument min_denominator is specified for functions analyzeVariants() or getSGVariantCounts(), estimates are set to NA unless at least one of x.S or x.E is equal or greater to the specified value.
Note that SGVariantCounts objects also store the raw count data. Count data can be used for statistical modeling, for example as suggested in section Testing for differential splice variant usage.
## N1 N2 N3 N4 T1 T2 T3 T4
## [1,] 0.88 0.56 0.6153846 0.5925926 0.93333333 0.90196078 0.8484848 0.95833333
## [2,] 0.12 0.44 0.3846154 0.4074074 0.06666667 0.09803922 0.1515152 0.04166667
Splice variants and estimates of relative usage can be visualized with function plotVariants.
## [updateObject] Validating the updated object ... OK
## [updateObject] Validating the updated object ... OK
plotVariants generates a two-panel figure similar to plotFeatures. The splice graph in the top panel illustrates the selected splice event. In this example, the event consists of two variants, corresponding to a skip or inclusion of the unannotated exon. The heatmap illustrates estimates of relative usage for each splice variant. Samples N2, N3 and N4 show evidence for transcripts that include the exon as well as transcripts that skip the exon. The remaining samples show little evidence for exon inclusion.
The functional consequences of a splice variant can be assessed by predicting its effect on protein-coding potential. Function predictVariantEffects() takes as input an SGVariants object with splice variants of interest, a set of annotated transcripts, and a matching reference genome provided as a BSgenome object.
## Warning in (function (seqlevels, genome, new_style) : cannot switch some hg19's
## seqlevels from UCSC to NCBI style
## variantID txName geneName RNA_change
## 1 2 uc002fjv.3 79791 r.88_89ins88+1798_88+1884
## 2 2 uc002fjw.3 79791 r.412_413ins412+1798_412+1884
## 3 2 uc010vot.2 79791 r.-105_-104ins-105+1798_-105+1884
## RNA_variant_type protein_change
## 1 insertion p.K29_L30insRINPRVKSGRFVKILPDYEHMAYRDVYTC
## 2 insertion p.K137_L138insRINPRVKSGRFVKILPDYEHMAYRDVYTC
## 3 insertion p.=
## protein_variant_type
## 1 in-frame_insertion
## 2 in-frame_insertion
## 3 no_change
The output is a data frame with each row describing the effect of a particular splice variant on an annotated protein-coding transcript. The effect of the variants is described following HGVS recommendations. In its current implementation, variant effect prediction is relatively slow and it is recommended to run predictVariantEffects() on select variants only.
Functions plotFeatures() and plotVariants() support many options for customizing figures. The splice graph in the top panel is plotted by function plotSpliceGraph, which can be called directly.
plotFeatures() includes multiple arguments for selecting features to be displayed. The following code block illustrates three different options for plotting the splice graph and expression levels for FBXO31 (Entrez ID 79791).
plotFeatures(sgfc_pred, geneID = 1)
plotFeatures(sgfc_pred, geneName = "79791")
plotFeatures(sgfc_pred, which = gr)
By default, the heatmap generated by plotFeatures() displays splice junctions. Alternatively, exon bins, or both exon bins and splice junctions can be displayed.
plotFeatures(sgfc_pred, geneID = 1, include = "junctions")
plotFeatures(sgfc_pred, geneID = 1, include = "exons")
plotFeatures(sgfc_pred, geneID = 1, include = "both")
Argument toscale controls which parts of the gene model are drawn to scale.
plotFeatures(sgfc_pred, geneID = 1, toscale = "gene")
plotFeatures(sgfc_pred, geneID = 1, toscale = "exon")
plotFeatures(sgfc_pred, geneID = 1, toscale = "none")
Heatmaps allow the visualization of expression values summarized for splice junctions and exon bins. Alternatively, per-base read coverages and splice junction counts can be visualized with function plotCoverage.
SGSeq does not implement statistical tests for differential splice variant usage. However, existing software packages such as DEXSeq (Anders, Reyes, and Huber 2012) and limma (Ritchie et al. 2015) can be used for this purpose. These packages allow infererence of differential exon usage within a gene (between groups of samples). In the SGSeq framework, this approach can be used to test for differential splice variant usage within a splice event, treating splice variants and splice events as exons and genes, respectively.
For these methods to be applicable, a single count is required for each splice variant. SGVariantCounts objects as described above store two counts for each splice variant, one for the 5′ and one for the 3′ end of the variant. These counts can be readily obtained from SGFeatureCounts objects, but are impractical for count-based differential testing. A single count for each variant (based on fragments compatible at either end of the variant) can be obtained from BAM files using function getSGVariantCounts(). The output is an SGVariantCounts object with additional assay countsVariant5pOr3p.
## class: SGVariantCounts
## dim: 2 8
## metadata(0):
## assays(6): countsVariant5p countsVariant3p ... countsVariant5pOr3p
## variantFreq
## rownames: NULL
## rowData names(20): from to ... variantType variantName
## colnames(8): N1 N2 ... T3 T4
## colData names(6): sample_name file_bam ... frag_length lib_size
Performing differential tests requires per-variant counts, unique identifiers for each variant, and a variable indicating how variants are grouped by events. All three can be obtained from the SGVariantCounts object.
Treating these three objects analogously to per-exon counts, exon identifiers and gene identifiers, they can be used to construct a DEXSeqDataSet object for use with DEXSeq or as input for function diffSplice() in combination with voom() for use with limma.
Functions analyzeFeatures() and analyzeVariants() wrap multiple analysis steps for convenience. Alternatively, the functions performing individual steps can be called directly. For example, the analysis based on de novo prediction can be performed as follows.
txf <- predictTxFeatures(si, gr)
sgf <- convertToSGFeatures(txf)
sgf <- annotate(sgf, txf_ucsc)
sgfc <- getSGFeatureCounts(si, sgf)
sgv <- findSGVariants(sgf)
sgvc <- getSGVariantCounts(sgv, sgfc)
predictTxFeatures() predicts features for each sample, merges features across samples, and performs filtering and processing of predicted terminal exons. predictTxFeatures() and getSGFeatureCounts() can also be run on individual samples (e.g. for distribution across a high-performance computing cluster). When using predictTxFeatures() for individual samples, with predictions intended to be merged later, run predictTxFeatures() with argument min_overhang = NULL to suppress processing of terminal exons. Then predictions can subsequently be merged and processed with functions mergeTxFeatures() and processTerminalExons(), respectively.
SGSeq supports parallelization across multiple cores on the same compute node. For functions that support parallelization, the number of cores can be specified with the cores argument. For most genome-wide analyses, memory requirements should not exceed 16 Gb per core. Memory issues encountered during feature prediction or read counting are likely due to individual genomic regions with high read coverage. By default, regions with high split-read complexity are skipped during prediction (controlled with argument max_complexity). For the counting step, memory issues may be resolved by excluding problematic regions (e.g. mitochondrial genes).
## 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] BSgenome.Hsapiens.UCSC.hg19_1.4.3
## [2] BSgenome_1.75.0
## [3] rtracklayer_1.66.0
## [4] BiocIO_1.17.0
## [5] TxDb.Hsapiens.UCSC.hg19.knownGene_3.2.2
## [6] GenomicFeatures_1.59.0
## [7] AnnotationDbi_1.69.0
## [8] SGSeq_1.41.0
## [9] SummarizedExperiment_1.36.0
## [10] Biobase_2.67.0
## [11] MatrixGenerics_1.19.0
## [12] matrixStats_1.4.1
## [13] Rsamtools_2.22.0
## [14] Biostrings_2.75.0
## [15] XVector_0.46.0
## [16] GenomicRanges_1.59.0
## [17] GenomeInfoDb_1.43.0
## [18] IRanges_2.41.0
## [19] S4Vectors_0.44.0
## [20] BiocGenerics_0.53.0
## [21] knitr_1.48
## [22] BiocStyle_2.35.0
##
## loaded via a namespace (and not attached):
## [1] KEGGREST_1.47.0 rjson_0.2.23 xfun_0.48
## [4] bslib_0.8.0 lattice_0.22-6 vctrs_0.6.5
## [7] tools_4.4.1 bitops_1.0-9 RUnit_0.4.33
## [10] curl_5.2.3 parallel_4.4.1 RSQLite_2.3.7
## [13] highr_0.11 blob_1.2.4 pkgconfig_2.0.3
## [16] Matrix_1.7-1 lifecycle_1.0.4 GenomeInfoDbData_1.2.13
## [19] compiler_4.4.1 codetools_0.2-20 htmltools_0.5.8.1
## [22] sys_3.4.3 buildtools_1.0.0 sass_0.4.9
## [25] RCurl_1.98-1.16 yaml_2.3.10 crayon_1.5.3
## [28] jquerylib_0.1.4 BiocParallel_1.41.0 DelayedArray_0.33.1
## [31] cachem_1.1.0 abind_1.4-8 digest_0.6.37
## [34] restfulr_0.0.15 maketools_1.3.1 fastmap_1.2.0
## [37] grid_4.4.1 cli_3.6.3 SparseArray_1.6.0
## [40] magrittr_2.0.3 S4Arrays_1.6.0 XML_3.99-0.17
## [43] UCSC.utils_1.2.0 bit64_4.5.2 rmarkdown_2.28
## [46] httr_1.4.7 igraph_2.1.1 bit_4.5.0
## [49] png_0.1-8 memoise_2.0.1 evaluate_1.0.1
## [52] rlang_1.1.4 glue_1.8.0 DBI_1.2.3
## [55] BiocManager_1.30.25 jsonlite_1.8.9 R6_2.5.1
## [58] GenomicAlignments_1.43.0 zlibbioc_1.52.0