Abstract
This vignette is a guide containing example code for performing real-life tasks. Importantly, it covers some functionality that were not covered in the Quick-Start vignette (because they are too computationally intensive to be reproducible in a vignette). Version 1.4.1
For instructions on installing and configuring SpliceWiz, please see the Quick-Start vignette.
library(SpliceWiz)
#> Loading required package: NxtIRFdata
#> SpliceWiz package loaded with 2 threads
#> Use setSWthreads() to set the number of SpliceWiz threads
First, define the path to the directory in which the reference should be stored. This directory will be made by SpliceWiz, but its parent directory must exist, otherwise an error will be returned.
Note that setting genome_path = "hg38"
will prompt SpliceWiz to use the default files for nonPolyA and Mappability exclusion references in the generation of its reference. Valid options for genome_path
are “hg38”, “hg19”, “mm10” and “mm9”.
buildRef(
reference_path = ref_path,
fasta = "genome.fa", gtf = "transcripts.gtf",
genome_type = "hg38"
)
buildRef()
first runs getResources()
, which prepares the genome and gene annotations by storing a compressed local copy in the resources
subdirectory of the given reference path. Specifically, a binary compressed version of the FASTA file (a.k.a. TwoBitFile), and a gzipped GTF file. If fasta
and/or gtf
are https or ftp links, the resources will be downloaded from the internet (which may take a while).
After local compressed versions of the genome and gene annotations are prepared, buildRef()
will proceed to generate the SpliceWiz reference.
Note that these two steps can be run separately. getResources()
will prepare local compressed copies of the FASTA / GTF resources without generating the SpliceWiz reference. Running buildRef()
, with reference_path
specifying where the resources were prepared previously with getResources()
, will perform the 2nd step (SpliceWiz reference generation) without needing to prepare the genome resources (in this case, set the parameters fasta = ""
and gtf = ""
).
As an example, the below steps:
getResources(
reference_path = ref_path,
fasta = "genome.fa",
gtf = "transcripts.gtf"
)
buildRef(
reference_path = ref_path,
fasta = "", gtf = "",
genome_type = "hg38"
)
is equivalent to this:
buildRef(
reference_path = ref_path,
fasta = "genome.fa",
gtf = "transcripts.gtf"
genome_type = "hg38"
)
To re-build and overwrite an existing reference, using the same resource annotations, set overwrite = TRUE
# Assuming hg38 genome:
buildRef(
reference_path = ref_path,
genome_type = "hg38",
overwrite = TRUE
)
If buildRef()
is run without setting overwrite = TRUE
, it will terminate if the file SpliceWiz.ref.gz
is found within the reference directory.
The following will first download the genome and gene annotation files from the online resource and store a local copy of it in a file cache, facilitated by BiocFileCache. Then, it uses the downloaded resource to create the SpliceWiz reference.
FTP <- "ftp://ftp.ensembl.org/pub/release-94/"
buildRef(
reference_path = ref_path,
fasta = paste0(FTP, "fasta/homo_sapiens/dna/",
"Homo_sapiens.GRCh38.dna.primary_assembly.fa.gz"),
gtf = paste0(FTP, "gtf/homo_sapiens/",
"Homo_sapiens.GRCh38.94.chr.gtf.gz"),
genome_type = "hg38"
)
AnnotationHub contains Ensembl references for many genomes. To browse what is available:
require(AnnotationHub)
#> Loading required package: AnnotationHub
#> Loading required package: BiocGenerics
#>
#> Attaching package: 'BiocGenerics'
#> 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, setdiff, sort, table,
#> tapply, union, unique, unsplit, which.max, which.min
#> Loading required package: BiocFileCache
#> Loading required package: dbplyr
ah <- AnnotationHub()
query(ah, "Ensembl")
#> AnnotationHub with 36085 records
#> # snapshotDate(): 2023-10-21
#> # $dataprovider: Ensembl, FANTOM5,DLRP,IUPHAR,HPRD,STRING,SWISSPROT,TREMBL,E...
#> # $species: Mus musculus, Homo sapiens, Sus scrofa, Rattus norvegicus, Danio...
#> # $rdataclass: TwoBitFile, GRanges, EnsDb, SQLiteFile, data.frame, OrgDb, li...
#> # additional mcols(): taxonomyid, genome, description,
#> # coordinate_1_based, maintainer, rdatadateadded, preparerclass, tags,
#> # rdatapath, sourceurl, sourcetype
#> # retrieve records with, e.g., 'object[["AH5046"]]'
#>
#> title
#> AH5046 | Ensembl Genes
#> AH5160 | Ensembl Genes
#> AH5311 | Ensembl Genes
#> AH5434 | Ensembl Genes
#> AH5435 | Ensembl EST Genes
#> ... ...
#> AH114088 | org.Mmu.eg.db.sqlite
#> AH114089 | org.Ce.eg.db.sqlite
#> AH114090 | org.Xl.eg.db.sqlite
#> AH114091 | org.Sc.sgd.db.sqlite
#> AH114092 | org.Dr.eg.db.sqlite
For a more specific query:
query(ah, c("Homo Sapiens", "release-94"))
#> AnnotationHub with 9 records
#> # snapshotDate(): 2023-10-21
#> # $dataprovider: Ensembl
#> # $species: Homo sapiens
#> # $rdataclass: TwoBitFile, GRanges
#> # additional mcols(): taxonomyid, genome, description,
#> # coordinate_1_based, maintainer, rdatadateadded, preparerclass, tags,
#> # rdatapath, sourceurl, sourcetype
#> # retrieve records with, e.g., 'object[["AH64628"]]'
#>
#> title
#> AH64628 | Homo_sapiens.GRCh38.94.abinitio.gtf
#> AH64629 | Homo_sapiens.GRCh38.94.chr.gtf
#> AH64630 | Homo_sapiens.GRCh38.94.chr_patch_hapl_scaff.gtf
#> AH64631 | Homo_sapiens.GRCh38.94.gtf
#> AH65744 | Homo_sapiens.GRCh38.cdna.all.2bit
#> AH65745 | Homo_sapiens.GRCh38.dna.primary_assembly.2bit
#> AH65746 | Homo_sapiens.GRCh38.dna_rm.primary_assembly.2bit
#> AH65747 | Homo_sapiens.GRCh38.dna_sm.primary_assembly.2bit
#> AH65748 | Homo_sapiens.GRCh38.ncrna.2bit
We wish to fetch “AH65745” and “AH64631” which contains the desired FASTA and GTF files, respectively. To build a reference using these resources:
Build-Reference-methods
will recognise the inputs of fasta
and gtf
as AnnotationHub resources if they begin with “AH”.
For human and mouse genomes, we highly recommend specifying genome_type
as the default mappability file is used to exclude intronic regions with repeat sequences from intron retention analysis. For other species, one could generate a SpliceWiz reference without this reference:
buildRef(
reference_path = ref_path,
fasta = "genome.fa", gtf = "transcripts.gtf",
genome_type = ""
)
If one wishes to prepare a Mappability Exclusion for species other than human or mouse, please see the Calculating Mappability Exclusions using STAR
section below.
For human and mouse genomes, gene ontology annotations are automatically generated. This is inferred by specifying genome_type
to the human or mouse genome. For other species, or to specify human/mouse, this should be specified in the ontologySpecies
parameter of buildRef()
.
Only Ensembl/orgDB resources are supported (for now). For a list of available species:
getAvailableGO()
#> [1] "Coffea arabica" "Anopheles gambiae"
#> [3] "Arabidopsis thaliana" "Bos taurus"
#> [5] "Canis familiaris" "Gallus gallus"
#> [7] "Pan troglodytes" "Escherichia coli"
#> [9] "Drosophila melanogaster" "Homo sapiens"
#> [11] "Mus musculus" "Sus scrofa"
#> [13] "Rattus norvegicus" "Macaca mulatta"
#> [15] "Caenorhabditis elegans" "Xenopus laevis"
#> [17] "Saccharomyces cerevisiae" "Danio rerio"
For example, to specify arabidopsis:
buildRef(
reference_path = ref_path,
fasta = "genome.fa", gtf = "transcripts.gtf",
genome_type = "",
ontologySpecies = "Arabidopsis thaliana"
)
To use STAR
to align FASTQ files, one must be using a system with STAR
installed. This software is not available in Windows. To check if STAR
is available:
ref_path = "./Reference"
# Ensure genome resources are prepared from genome FASTA and GTF file:
if(!dir.exists(file.path(ref_path, "resource"))) {
getResources(
reference_path = ref_path,
fasta = "genome.fa",
gtf = "transcripts.gtf"
)
}
# Generate a STAR genome reference:
STAR_BuildRef(
reference_path = ref_path,
n_threads = 8
)
Note that, by default, STAR_BuildRef
will store the STAR genome reference in the STAR
subdirectory within reference_path
. To override this setting, set the STAR_ref_path
parameter to a directory path of your choice, e.g.:
STAR_BuildRef(
reference_path = ref_path,
STAR_ref_path = "/path/to/another/directory",
n_threads = 8
)
Sometimes, one might wish to build a genome annotation without first specifying the gene annotations. Reasons one might want to do this include:
We can use STAR_buildGenome
to do this:
# Generate a STAR genome reference:
STAR_buildGenome(
reference_path = ref_path,
STAR_ref_path = "/path/to/hg38"
n_threads = 8
)
This STAR reference is derived from the genome FASTA file but not the gene annotation GTF file. Prior to alignment, additional parameters need to be supplied (which should take 5 minutes). These include:
reference_path
parameterTo generate an on-the-fly (i.e., alignment-ready) STAR reference from a genome-derived reference:
STAR_new_ref <- STAR_loadGenomeGTF(
reference_path = ref_path,
STAR_ref_path = "/path/to/hg38",
STARgenome_output = file.path(tempdir(), "STAR"),
n_threads = 8,
sjdbOverhang = 100,
extraFASTA = "./ercc.fasta"
)
The path to the on-the-fly reference is specified by the return value (STAR_new_ref
in the above example).
As already explained, this step allows a single STAR reference to be built for each species, which can be adapted for different projects based on their specific technical specifications (e.g. different read length can be adapted by setting different sjdbOverhang
, or any spike-ins by setting the spike-in FASTA using extraFASTA
).
Genomes contain regions of low mappability (i.e. areas which are difficult for reads or fragments to align to). A common computational cause of low mappability include repeat sequences. IRFinder uses an empirical method to determine regions of low mappability, which we adopted in SpliceWiz. These resources are used automatically when generating the SpliceWiz reference and setting the genome_type
to supported genomes (hg38, hg19, mm10, mm9). For other species, one may wish to generate their own annotations of low mappability regions using the STAR aligner.
The STAR_mappability
wrapper function will use the STAR aligner to calculate regions of low mappability within the given genome.
STAR_mappability(
reference_path = ref_path,
STAR_ref_path = file.path(ref_path, "STAR"),
map_depth_threshold = 4,
n_threads = 8,
read_len = 70,
read_stride = 10,
error_pos = 35
)
In the above example, STAR_mappability()
will use the given STAR reference (inside the STAR_ref_path
directory), and the genome found within the reference_path
SpliceWiz reference, to generate synthetic reads.
read_len
specifies the length of these synthetic reads (default 70
)read_stride
specifies the nucleotide distance between adjacent synthetic reads (default 10
). These will be generated with alternate +
/ -
stranderror_pos
introduces a single nucleotide error at the specified position (default 35
), which will generate an SNP at the center of the 70-nt synthetic read.These synthetic reads will then be aligned back to the STAR genome to create a BAM file, which is later processed to measure the coverage depth of the genome by these synthetic reads.
Finally, regions with coverage depth of map_depth_threshold
or below will be defined as regions of “low mappability”. In the above example, 70-nt reads of 10-nt stride will produce synthetic reads such that each nucleotide is expected to have a coverage of 70 / 10 = 7
nucleotides. A coverage of 4
nucleotides or less equates to a coverage of < ~60% of expected depth.
If STAR
is available on the same computer or server where R/RStudio is being run, we can use the one-line function buildFullRef
. This function will:
getResources
)STAR_BuildRef
)STAR_mappability
)buildRef
)This step is recommended when one wishes to build a non-human/mouse genome in a single step, including generating low-mappability regions to exclude measuring IR events with low mappability.
buildFullRef(
reference_path = ref_path,
fasta = "genome.fa", gtf = "transcripts.gtf",
genome_type = "",
use_STAR_mappability = TRUE,
n_threads = 8
)
n_threads
specify how many threads should be used to build the STAR reference and to calculate the low mappability regions
If STAR
is not available, Rsubread
is available on Bioconductor for alignment and can be used to perform mappability calculations. The example code in the manual is displayed here for convenience, to demonstrate how this would be done:
require(Rsubread)
# (1a) Creates genome resource files
ref_path <- file.path(tempdir(), "Reference")
getResources(
reference_path = ref_path,
fasta = chrZ_genome(),
gtf = chrZ_gtf()
)
# (1b) Systematically generate reads based on the SpliceWiz example genome:
generateSyntheticReads(
reference_path = ref_path
)
# (2) Align the generated reads using Rsubread:
# (2a) Build the Rsubread genome index:
subreadIndexPath <- file.path(ref_path, "Rsubread")
if(!dir.exists(subreadIndexPath)) dir.create(subreadIndexPath)
Rsubread::buildindex(
basename = file.path(subreadIndexPath, "reference_index"),
reference = chrZ_genome()
)
# (2b) Align the synthetic reads using Rsubread::subjunc()
Rsubread::subjunc(
index = file.path(subreadIndexPath, "reference_index"),
readfile1 = file.path(ref_path, "Mappability", "Reads.fa"),
output_file = file.path(ref_path, "Mappability", "AlignedReads.bam"),
useAnnotation = TRUE,
annot.ext = chrZ_gtf(),
isGTF = TRUE
)
# (3) Analyse the aligned reads in the BAM file for low-mappability regions:
calculateMappability(
reference_path = ref_path,
aligned_bam = file.path(ref_path, "Mappability", "AlignedReads.bam")
)
# (4) Build the SpliceWiz reference using the calculated Mappability Exclusions
buildRef(ref_path)
Note that the default output file for calculateMappability()
(step 3) is Mappability/MappabilityExclusion.bed.gz
found within the reference_path
directory. Then buildRef()
(step 4) will automatically use this file, regardless of the genome_type
parameter. The exception is if MappabilityRef
parameter is set to a different file.
This conveniences users to generate their own human/mouse mappability files but use the default non-polyA reference, e.g.:
First, remember to check that STAR is available via command line:
STAR_alignReads(
fastq_1 = "sample1_1.fastq", fastq_2 = "sample1_2.fastq",
STAR_ref_path = file.path(ref_path, "STAR"),
BAM_output_path = "./bams/sample1",
n_threads = 8,
trim_adaptor = "AGATCGGAAG"
)
Note that by default, STAR_alignReads()
will “trim” Illumina adapters (in fact they will be soft-clipped using STAR’s --clip3pAdapterSeq
option). To disable this feature, set trim_adapter = ""
in the STAR_alignReads()
function.
Experiment <- data.frame(
sample = c("sample_A", "sample_B"),
forward = file.path("raw_data", c("sample_A", "sample_B"),
c("sample_A_1.fastq", "sample_B_1.fastq")),
reverse = file.path("raw_data", c("sample_A", "sample_B"),
c("sample_A_2.fastq", "sample_B_2.fastq"))
)
STAR_alignExperiment(
Experiment = Experiment,
STAR_ref_path = file.path("Reference_FTP", "STAR"),
BAM_output_path = "./bams",
n_threads = 8,
two_pass = FALSE
)
To use two-pass mapping, set two_pass = TRUE
. We recommend disabling this feature, as one-pass mapping is adequate in typical-use cases. Two-pass mapping is recommended if one expects a large number of novel splicing events or if the gene annotations (of transcript isoforms) is likely to be incomplete. Additionally, two-pass mapping is highly memory intensive and should be reserved for systems with high memory resources.
SpliceWiz can identify sequencing FASTQ files recursively from a given directory. It assumes that forward and reverse reads are suffixed as _1
and _2
, respectively. Users can choose to identify such files using a specified file extension. For example, to recursively identify FASTQ files of the format {sample}_1.fq.gz
and {sample}_2.fq.gz
, use the following:
# Assuming sequencing files are named by their respective sample names
fastq_files <- findFASTQ(
sample_path = "./sequencing_files",
paired = TRUE,
fastq_suffix = ".fq.gz", level = 0
)
For gzipped fastq files, fastq_suffix
should be ".fq.gz"
or ".fastq.gz"
. For uncompressed fastq files, it should be ".fq"
or ".fastq"
. Please check your files in order to correctly set this option.
findFASTQ()
will return a 2- or 3-column data frame (depending if paired
was set to FALSE
or TRUE
, respectively). The first column is the sample name (the file name, if level = 0
, or the parent directory name, if level = 1
). The subsequent columns are the paths of the forward and reverse reads.
The data.frame returned by the findFASTQ()
function can be parsed into the STAR_alignExperiment
function. This will align all samples contained in the data.frame parsed via the Experiment
parameter.
STAR_alignExperiment(
Experiment = fastq_files,
STAR_ref_path = file.path("Reference_FTP", "STAR"),
BAM_output_path = "./bams",
n_threads = 8,
two_pass = FALSE
)
Note that, if a directory contains multiple forward and reverse FASTQ files, they will be aligned to the same BAM file. This can be done by setting level = 1
in the findFASTQ()
function, resulting in multiple rows with the same sample name.
To conveniently find all BAM files recursively in a given path:
This convenience function returns the putative sample names, either from BAM file names themselves (level = 0
), or from the names of their parent directories (level = 1
).
First, ensure that a SpliceWiz reference has been generated using the buildRef()
function. This reference should be parsed into the reference_path
parameter of the processBAM()
function.
To run processBAM()
using 4 OpenMP threads:
# assume SpliceWiz reference has been generated in `ref_path` using the
# `buildRef()` function.
processBAM(
bamfiles = bams$path,
sample_names = bams$sample,
reference_path = ref_path,
output_path = "./pb_output",
n_threads = 4,
useOpenMP = TRUE
)
Sometimes one may wish to create a COV file from a BAM file without running processBAM()
. One reason might be because a SpliceWiz reference is not available.
To convert a list of BAM files, run BAM2COV()
. This is a function structurally similar to processBAM()
but without the need to give the path to the SpliceWiz reference:
BAM2COV(
bamfiles = bams$path,
sample_names = bams$sample,
output_path = "./cov_output",
n_threads = 4,
useOpenMP = TRUE
)
Sometimes, users may wish to convert COV files to BigWig. One common reason may be to generate strand-specific coverage to compare with BigWig files on IGV.
For example, to generate a BigWig file containing reads on the negative strand:
se <- SpliceWiz_example_NxtSE()
cov_file <- covfile(se)[1]
cov_negstrand <- getCoverage(cov_file, strand = "-")
bw_file <- file.path(tempdir(), "sample_negstrand.bw")
rtracklayer::export(cov_negstrand, bw_file, "bw")
SpliceWiz processes BAM files using OpenMP-based parallelisation (multi-threading), using our ompBAM C++ library (available via the ompBAM Bioconductor package). The advantage of using this approach (instead of processing multiple BAM files each using a single thread) is that the latter approach uses a lot more memory. Our OpenMP-based approach processes BAM files one at a time, avoiding the memory cost when analysing multiple BAM files simultaneously.
Note that, by default, processBAM
and BAM2COV
will use OpenMP where available (which is natively supported on Windows and Linux). For MacOS, if OpenMP is not available, these functions will use BiocParallel’s MulticoreParam
to multi-thread process BAM files (1 BAM per thread). Beware that this may take a lot of RAM! (Typically 5-10 Gb per BAM file). We highly suggest considering installing OpenMP libraries on MacOS, as this will lower RAM usage.
Assuming the SpliceWiz reference is in ref_path
, after running processBAM()
as shown in the previous section, use the convenience function findSpliceWizOutput()
to tabulate a list of samples and their corresponding processBAM()
outputs:
This data.frame can be directly used to run collateData
:
novelSplicing = TRUE
. See the Quick-Start vignette for more details about the various parameters associated with novel splicing detection.collateData(
Experiment = expr,
reference_path = ref_path,
output_path = "./NxtSE_output_novelSplicing",
novelSplicing = TRUE
)
Then, the collated data can be imported as a NxtSE
object, which is an object that inherits SummarizedExperiment
and has specialized containers to hold additional data required by SpliceWiz.
Please refer to SpliceWiz: Quick-Start vignette for worked examples using the example dataset.
sessionInfo()
#> R version 4.3.2 Patched (2023-11-13 r85521)
#> Platform: x86_64-pc-linux-gnu (64-bit)
#> Running under: Ubuntu 22.04.3 LTS
#>
#> Matrix products: default
#> BLAS: /home/biocbuild/bbs-3.18-bioc/R/lib/libRblas.so
#> LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.10.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] stats graphics grDevices utils datasets methods base
#>
#> other attached packages:
#> [1] AnnotationHub_3.10.0 BiocFileCache_2.10.1 dbplyr_2.4.0
#> [4] BiocGenerics_0.48.1 SpliceWiz_1.4.1 NxtIRFdata_1.8.0
#>
#> loaded via a namespace (and not attached):
#> [1] RColorBrewer_1.1-3 jsonlite_1.8.8
#> [3] magrittr_2.0.3 rmarkdown_2.25
#> [5] fs_1.6.3 BiocIO_1.12.0
#> [7] zlibbioc_1.48.0 vctrs_0.6.5
#> [9] memoise_2.0.1 Rsamtools_2.18.0
#> [11] DelayedMatrixStats_1.24.0 RCurl_1.98-1.13
#> [13] webshot_0.5.5 progress_1.2.3
#> [15] htmltools_0.5.7 S4Arrays_1.2.0
#> [17] curl_5.2.0 Rhdf5lib_1.24.1
#> [19] SparseArray_1.2.3 rhdf5_2.46.1
#> [21] shinyFiles_0.9.3 sass_0.4.8
#> [23] bslib_0.6.1 htmlwidgets_1.6.4
#> [25] plotly_4.10.3 cachem_1.0.8
#> [27] GenomicAlignments_1.38.0 iterators_1.0.14
#> [29] mime_0.12 lifecycle_1.0.4
#> [31] pkgconfig_2.0.3 Matrix_1.6-4
#> [33] R6_2.5.1 fastmap_1.1.1
#> [35] GenomeInfoDbData_1.2.11 MatrixGenerics_1.14.0
#> [37] shiny_1.8.0 digest_0.6.33
#> [39] colorspace_2.1-0 patchwork_1.1.3
#> [41] AnnotationDbi_1.64.1 S4Vectors_0.40.2
#> [43] GenomicRanges_1.54.1 RSQLite_2.3.4
#> [45] seriation_1.5.4 filelock_1.0.3
#> [47] fansi_1.0.6 httr_1.4.7
#> [49] abind_1.4-5 compiler_4.3.2
#> [51] withr_2.5.2 bit64_4.0.5
#> [53] BiocParallel_1.36.0 viridis_0.6.4
#> [55] DBI_1.2.0 heatmaply_1.5.0
#> [57] dendextend_1.17.1 HDF5Array_1.30.0
#> [59] R.utils_2.12.3 rappdirs_0.3.3
#> [61] DelayedArray_0.28.0 rjson_0.2.21
#> [63] tools_4.3.2 interactiveDisplayBase_1.40.0
#> [65] httpuv_1.6.13 fst_0.9.8
#> [67] R.oo_1.25.0 glue_1.6.2
#> [69] restfulr_0.0.15 rhdf5filters_1.14.1
#> [71] promises_1.2.1 grid_4.3.2
#> [73] generics_0.1.3 gtable_0.3.4
#> [75] BSgenome_1.70.1 ca_0.71.1
#> [77] R.methodsS3_1.8.2 tidyr_1.3.0
#> [79] hms_1.1.3 data.table_1.14.10
#> [81] xml2_1.3.6 utf8_1.2.4
#> [83] XVector_0.42.0 foreach_1.5.2
#> [85] BiocVersion_3.18.1 pillar_1.9.0
#> [87] genefilter_1.84.0 later_1.3.2
#> [89] splines_4.3.2 dplyr_1.1.4
#> [91] lattice_0.22-5 ompBAM_1.6.0
#> [93] survival_3.5-7 rtracklayer_1.62.0
#> [95] bit_4.0.5 annotate_1.80.0
#> [97] tidyselect_1.2.0 registry_0.5-1
#> [99] Biostrings_2.70.1 knitr_1.45
#> [101] rhandsontable_0.3.8 gridExtra_2.3
#> [103] IRanges_2.36.0 SummarizedExperiment_1.32.0
#> [105] stats4_4.3.2 xfun_0.41
#> [107] shinydashboard_0.7.2 Biobase_2.62.0
#> [109] matrixStats_1.2.0 pheatmap_1.0.12
#> [111] DT_0.31 stringi_1.8.3
#> [113] lazyeval_0.2.2 yaml_2.3.8
#> [115] shinyWidgets_0.8.0 evaluate_0.23
#> [117] codetools_0.2-19 tibble_3.2.1
#> [119] BiocManager_1.30.22 cli_3.6.2
#> [121] xtable_1.8-4 munsell_0.5.0
#> [123] jquerylib_0.1.4 Rcpp_1.0.11
#> [125] GenomeInfoDb_1.38.5 png_0.1-8
#> [127] XML_3.99-0.16 parallel_4.3.2
#> [129] fstcore_0.9.18 ellipsis_0.3.2
#> [131] ggplot2_3.4.4 assertthat_0.2.1
#> [133] blob_1.2.4 prettyunits_1.2.0
#> [135] sparseMatrixStats_1.14.0 bitops_1.0-7
#> [137] viridisLite_0.4.2 scales_1.3.0
#> [139] purrr_1.0.2 crayon_1.5.2
#> [141] rlang_1.1.2 rvest_1.0.3
#> [143] TSP_1.2-4 KEGGREST_1.42.0