CRISPR-Cas nucleases and their derivatives, such as nickases and gene expression regulators, have become the most popular tools for genome modification and gene expression regulation in both model organisms and human cells (Mali et al. 2013; Hsu et al. 2013). The CRISPR-Cas system uses a guide RNA (gRNA) to direct the Cas nuclease to target DNA sequence that is adjacent to a Protospace Adjacent Motif (PAM). The Cas enzyme then creates a double-strand break at the target site, initiating the gene modification. Following the cut, various editing outcomes, such as short insertions and deletions (INDELs) as well as point mutations, can occur through different DNA damage repair mechanisms (Chen et al. 2019). In the most widely used CRISPR system from Streptococcus pyogenes, the gRNA consists of 20 nucleotides, with the preferred PAM being “NGG” (or “NAG” with reduced activity).
A key limitation of CRISPR systems is their potential to cleave DNA sequences that do no perfectly match the gRNA target, known as “off-target effects”. Therefore, designing gRNAs with low off-target activity is crucial for the effective application of CRISPR systems (Lihua Julie Zhu 2015). Several strategies can help minimize off-target effects:
We developed the CRISPRseek package to assist with all the above steps. Key features include:
The CRISPRseek package generates several reports. Key output includes:
The CRISPRseek package is highly flexible, allowing users to customize gRNA and PAM sequence requirements for different CRISPR systems across various bacterial species. It can also be easily extended to incorporate improved weight matrices or scoring methods for off-target analysis as new experimental and computational results emerge.
Key functions implemented in the CRISPRseek,
along with their roles, input parameters, and outputs, are illustrated
in the diagram below. For convenience, we highly recommend using the
one-stop wrapper function, offtargetAnalysis()
, which
streamlines the entire gRNA identification and off-target analysis
workflow.
Analytical workflow and core functions in CRISPRseek
In this section, we will demonstrate different gRNA design scenarios
using the offtargetAnalysis()
function. Some of its core
parameters are described below. Type ?offtargetAnalysis
for
detailed description of all supported parameters.
DNAStringSet
object
containing sequences to be searched for potential gRNAs.inputFilePath
. Set to FALSE if the input file already
contains user-selected gRNAs plus PAM.BSgenome
object containing the target genome
sequence, used for off-target search.BSgenomeName
. Specifies the path to a
custom target genome file in FASTA format, used for off-target search.
It is applicable when BSgenomeName
is NOT set. When
genomeSeqFile
is set, the annotateExon
,
txdb
, and orgAnn
parameters will be
ignored.TxDb
object containing organism-specific annotation
data, required for annotateExon
.OrgDb
object containing organism-specific annotation
mapping information, required for annotateExon
.To annotate the resulting off-targets to nearby exons, the following parameters are required: “annotateExon”, “BSgenomeName”, “txdb”, and “orgAnn”.
available.genomes()
function in the BSgenome
package. Common examples include BSgenome.Hsapiens.UCSC.hg19
(hg19), BSgenome.Hsapiens.UCSC.hg38
(hg38), BSgenome.Mmusculus.UCSC.mm10
(mm10), BSgenome.Celegans.UCSC.ce6
(ce6).TxDb
objects, search for
annotation packages starting with “TxDb” at Bioconductor.
Examples include TxDb.Hsapiens.UCSC.hg19.knownGene
(for hg19) and TxDb.Mmusculus.UCSC.mm10.knownGene
(for mm10). To build custom TxDb
, refer to the GenomicFeatures
and txdbmaker
packages.OrgDb
packages, search for
“OrgDb” at Bioconductor.
Examples include org.Hs.eg.db
(for human) and org.Mm.eg.db
(for mouse).The offTargetAnalysis()
function offers two input
options:
By default, the parameter findgRNAs = TRUE
directs the
function to accept a sequence file (via inputFilePath
),
search for potential gRNAs, and then perform off-target analysis.
Alternatively, if you already have a list of user-designed gRNAs, you
can provide them via inputFilePath
as well, set
findgRNAs = FALSE
, and the function will perform off-target
analysis without searching for gRNAs.
The following examples use a raw sequence as input and, therefore,
set findgRNAs = TRUE
.
By default, the offTargetAnalysis()
function will
identify all potential gRNAs in the given input sequence, perform
off-target analysis, and annotate for all identified off-targets. This
generates the most comprehensive reports, but comes with the trade-off
of the slowest running time. Additionally, we need to load the necessary
annotation packages first:
library(CRISPRseek)
library(BSgenome.Hsapiens.UCSC.hg38)
library(TxDb.Hsapiens.UCSC.hg38.knownGene)
library(org.Hs.eg.db)
Note that we set chromToSearch = "chrX"
to restrict the
off-target search to the “chrX” to speed up the analysis.
inputFilePath <- system.file('extdata', 'inputseq.fa', package = 'CRISPRseek')
outputDir <- getwd()
res <- offTargetAnalysis(inputFilePath = inputFilePath,
BSgenomeName = Hsapiens,
chromToSearch= "chrX",
# Annotation packages required for annotation.
txdb = TxDb.Hsapiens.UCSC.hg38.knownGene,
orgAnn = org.Hs.egSYMBOL,
outputDir = outputDir,
overwrite = TRUE)
head(res$summary[, c("names", "gRNAsPlusPAM", "gRNAefficacy")])
## names gRNAsPlusPAM gRNAefficacy
## 1 Hsap_GATA1_ex2_gR17f CCAGTTTGTGGATCCTGCTCNGG 0.10654794
## 2 Hsap_GATA1_ex2_gR20r GGTGTGGAGGACACCAGAGCNGG 0.14209108
## 3 Hsap_GATA1_ex2_gR39f GTGTCCTCCACACCAGAATCNGG 0.06962952
## 4 Hsap_GATA1_ex2_gR40f TGTCCTCCACACCAGAATCANGG 0.26188176
## 5 Hsap_GATA1_ex2_gR7r CCAGAGCAGGATCCACAAACNGG 0.10171983
## name OffTargetSequence alignment score
## 81 Hsap_GATA1_ex2_gR40f TGTCCTTGACACAAGAATCATGT ......TG....A....... 0.7
## 171 Hsap_GATA1_ex2_gR20r GGTCTGGAGGACAGCACAGCCTG ...C.........G..C... 0.2
## 191 Hsap_GATA1_ex2_gR20r TGTGTGTAGGACTCCAGAGCAAG T.....T.....T....... 0.8
## 16 Hsap_GATA1_ex2_gR20r GGTGTGGAGGAGATGAGAGCATG ...........G.TG..... 0.0
## 23 Hsap_GATA1_ex2_gR7r CCAGAGCAGGATCCACAAACTGG .................... 100.0
## 9 Hsap_GATA1_ex2_gR40f TGTCCTCCACACCAGAATCAGGG .................... 100.0
## gRNAefficacy
## 81 0.2068879
## 171 0.2487753
## 191 0.1297587
## 16 0.1521223
## 23 0.1017198
## 9 0.2618818
To skip the annotation step, simply set
annotateExon = FALSE
:
For a quicker gRNA search, you can call the function with
chromToSearch = ""
, which disables the search of identified
gRNAs against the reference genome. This will significantly reduce the
running time.
res <- offTargetAnalysis(inputFilePath = inputFilePath,
BSgenomeName = Hsapiens, # optional
chromToSearch = "",
outputDir = outputDir,
overwrite = TRUE)
res
## DNAStringSet object of length 5:
## width seq names
## [1] 23 CCAGTTTGTGGATCCTGCTCTGG Hsap_GATA1_ex2_gR17f
## [2] 23 GTGTCCTCCACACCAGAATCAGG Hsap_GATA1_ex2_gR39f
## [3] 23 TGTCCTCCACACCAGAATCAGGG Hsap_GATA1_ex2_gR40f
## [4] 23 GGTGTGGAGGACACCAGAGCAGG Hsap_GATA1_ex2_gR20r
## [5] 23 CCAGAGCAGGATCCACAAACTGG Hsap_GATA1_ex2_gR7r
Note that setting chromToSearch = ""
also disables the
calculation of gRNA efficacy, which measures how effectively a gRNA
facilitates cleavage at the target site. To enable the calculation of
gRNA efficacy at the on-target site, specify the chromosome where the
on-targets are located (e.g., chromToSearch = "chrX"
) and
set max.mismatch = 0
to ensure no off-target analysis is
performed.
res <- offTargetAnalysis(inputFilePath = inputFilePath,
BSgenomeName = Hsapiens,
annotateExon = FALSE,
chromToSearch = "chrX",
max.mismatch = 0,
outputDir = outputDir,
overwrite = TRUE)
head(res$summary[, c("names", "gRNAsPlusPAM", "gRNAefficacy")])
## names gRNAsPlusPAM gRNAefficacy
## 1 Hsap_GATA1_ex2_gR17f CCAGTTTGTGGATCCTGCTCNGG 0.10654794
## 2 Hsap_GATA1_ex2_gR20r GGTGTGGAGGACACCAGAGCNGG 0.14209108
## 3 Hsap_GATA1_ex2_gR39f GTGTCCTCCACACCAGAATCNGG 0.06962952
## 4 Hsap_GATA1_ex2_gR40f TGTCCTCCACACCAGAATCANGG 0.26188176
## 5 Hsap_GATA1_ex2_gR7r CCAGAGCAGGATCCACAAACNGG 0.10171983
Four rule sets are currently supported for calculating gRNA efficacy
(rule.set = c("Root_RuleSet1_2014", "Root_RuleSet2_2016", "CRISPRscan", "DeepCpf1")
)
(Doench et al.
2014, 2016; Moreno-Mateos et al. 2015; Kim et al. 2018). By default,
“Root_RuleSet_2014” is be used. To use “Root_RuleSet2_2016”, first
install python 2.7 via Anaconda (As RuleSet2 is implemented with python
2.7), then install the following Python packages: scikit-learn 0.16.1,
pickle, pandas, numpy, and scipy. Afterward, in an R session, run the
following commands:
Sys.setenv(PATH = paste("/anaconda2/bin", Sys.getenv("PATH"), sep = ":"))
system("python --version") # Should output Python 2.7.15.
If rule.set = "CRISPRscan"
, must also specify the
corresponding parameters for baseBefroegRNA
,
baseAfterPAM
, and featureWeightMatrixFile
to
ensure correct efficacy calculation:
m <- system.file("extdata", "Morenos-Mateo.csv", package = "CRISPRseek")
res <- offTargetAnalysis(inputFilePath = inputFilePath,
BSgenomeName = Hsapiens,
annotateExon = FALSE,
chromToSearch = "chrX",
max.mismatch = 0,
rule.set = "CRISPRscan",
baseBeforegRNA = 6,
baseAfterPAM = 6,
featureWeightMatrixFile = m,
outputDir = outputDir,
overwrite = TRUE)
If you would like to search for off-targets in custom reference
genomes (rather than using a BSgenome
), you can specify the
path to the custom reference sequence file via the
genomeSeqFile
argument. Please note that, when a custom
reference is used, arguments annotateExon
,
BSgenomeName
, txdb
, and
fetchSequence
will be ignored.
inputFilePath2 <- system.file("extdata", "inputseqWithoutBSgenome.fa", package = "CRISPRseek")
genomeSeqFile <- system.file("extdata", "genomeSeq.fasta", package = "CRISPRseek")
res <- offTargetAnalysis(inputFilePath = inputFilePath2,
genomeSeqFile = genomeSeqFile,
outputDir = outputDir,
overwrite = TRUE)
head(res$summary)
The CRISPRseek
supports searching for off-targets with bulges (both RNA bulges and DNA
bulges) by integrating the Cas-OFFinder (Bae, Park, and Kim 2014) tool. To enable
this, you can call the master function offTargetAnalysis()
with the argument findOffTargetsWithBulge = TRUE
.
There are three parameters specific to bulge searches, with the following default values:
method.findOffTargetsWithBulge = "CasOFFinder_v3.0.0b3"
DNA_bulge = 2
RNA_bulge = 2
Note that, currently, method.findOffTargetsWithBulge
only supports “CasOFFinder_v3.0.0b3”, which generates results that may
differ from those produced by “CasOFFinder_v2”, For more details, refer
to this link. To
use “CasOFFinder” on Linux or Windows, you may need to install
“pocl-opencl-icd” if it is not already installed.
if (interactive()) {
res <- offTargetAnalysis(inputFilePath = inputFilePath,
findOffTargetsWithBulge = TRUE,
BSgenomeName = Hsapiens,
chromToSearch = "chrX",
annotateExon = FALSE,
outputDir = outputDir,
overwrite = TRUE)
head(res$offtarget[, c("name", "gRNAPlusPAM_bulge", "OffTargetSequence_bulge", "n.RNABulge", "n.DNABulge", "alignment")])
}
The columns relevant to the bulge search are described below. If no bulge is detected, both the “gRNAPlusPAM_bulge” and “OffTargetSequence_bulge” columns for that row will both be empty, while the “n.RNABulge” and “n.DNABulge” values will both be 0.
Alternatively, you can use the wrapper function
getOfftargetWithBulge()
to directly output the Cas-OFFinder
results. The function accepts input in the form of either a
DNAStringSet
object from findgRNA()
, or a
list
object from offTargetAnalysis()
.
Note that getOfftargetWithBulge()
currently supports two
versions of Cas-OFFinder: “2.4.1” and “3.0.0b3”, specified through the
cas_offinder_version
argument, with the default set to
“2.4.1”. Type ?getOfftargetWithBulges
for more
examples.
The scoring.method
argument in
offTargetAnalysis()
determines how off-target scores are
calculated. By default, scoring.method = "Hsu-Zhuang"
,
which models the effect of mismatch position on cutting frequency. To
account for both mismatch position and type, you can switch to using the
CFD score (Doench et
al. 2016). To enable this, set
scoring.method = "CFDscore"
and
PAM.pattern = "NNG$|NGN$"
.
res <- offTargetAnalysis(inputFilePath = inputFilePath,
BSgenomeName = Hsapiens,
chromToSearch = "chrX",
annotateExon = FALSE,
scoring.method = "CFDscore",
PAM.pattern = "NNG$|NGN$",
outputDir = outputDir,
overwrite = TRUE)
head(res$offtarget[, c("name", "gRNAPlusPAM", "score", "alignment")])
## name gRNAPlusPAM score alignment
## 23 Hsap_GATA1_ex2_gR7r CCAGAGCAGGATCCACAAACNGG 1.000000 ....................
## 22 Hsap_GATA1_ex2_gR7r CCAGAGCAGGATCCACAAACNGG 0.032997 ....T.......T....T..
## 9 Hsap_GATA1_ex2_gR40f TGTCCTCCACACCAGAATCANGG 1.000000 ....................
## 14 Hsap_GATA1_ex2_gR40f TGTCCTCCACACCAGAATCANGG 0.022774 G...A...T...........
## 10 Hsap_GATA1_ex2_gR40f TGTCCTCCACACCAGAATCANGG 0.022096 ....T......AT.......
## 13 Hsap_GATA1_ex2_gR40f TGTCCTCCACACCAGAATCANGG 0.021065 .......A..........G.
You can filter the output to include only the desired gRNAs, such as paired gRNAs or those with specific restriction enzyme recognition sites, by adjusting the following parameters.
findPairedgRNAOnly = FALSE
:
min.gap
(defaults to 0) and
max.gap
(defaults to 20), inclusive. And the reverse gRNA
must be positioned before the forward gRNA.findgRNAsWithREcutOnly = FALSE
:
REpatternFile
parameter to specify the file
containing the restriction enzyme recognition patterns.Both parameters can be set to TRUE simultaneously to report only paired gRNAs, where at least one gRNA in the pair overlaps with a restriction enzyme recognition site.
res <- offTargetAnalysis(inputFilePath = inputFilePath,
BSgenomeName = Hsapiens,
annotateExon = FALSE,
chromToSearch = "chrX",
findPairedgRNAOnly = TRUE,
findgRNAsWithREcutOnly = TRUE,
outputDir = outputDir,
overwrite = TRUE)
head(res$summary[, c("names", "gRNAsPlusPAM", "gRNAefficacy", "PairedgRNAName", "REname")])
## names gRNAsPlusPAM gRNAefficacy
## 1 Hsap_GATA1_ex2_gR39f GTGTCCTCCACACCAGAATCNGG 0.06962952
## 2 Hsap_GATA1_ex2_gR40f TGTCCTCCACACCAGAATCANGG 0.26188176
## 3 Hsap_GATA1_ex2_gR7r CCAGAGCAGGATCCACAAACNGG 0.10171983
## PairedgRNAName REname
## 1 Hsap_GATA1_ex2_gR7r BslI HinfI TfiI
## 2 Hsap_GATA1_ex2_gR7r BslI HinfI TfiI
## 3 Hsap_GATA1_ex2_gR39f Hsap_GATA1_ex2_gR40f BslI PflMI
Searching for gRNAs in long input sequences (> 200 kb) may be
slow. To improve performance, set annotatePaired = FALSE
,
and enable multicore processing by setting
enable.multicore = TRUE
and adjusting
n.cores.max
. Additionally, we recommend splitting very long
sequences into smaller chunks and analyzing each sub-sequence
separately. Special thanks to Alex Williams for sharing this use case.
Finally, please ensure that repeat-masked sequences are used as
input.
To identify gRNAs that preferentially target one allele, the function
compare2Sequences()
, which takes two input files, can be
used. In the following example, file “rs362331C.fa” and “rs362331T.fa”
represent two input files that contain sequences differing by a single
nucleotide polymorphism (SNP). The results are saved in the file
“scoresFor2InputSequences.xlsx”. The output file lists all possible gRNA
sequences for each of the two input files and provides a cleavage score
for each of the two input sequences. To preferentially target one
allele, select gRNA sequences that have the lowest score for the other
allele. Selected gRNAs can then by examined for off-targets using
offTargetAnalysis()
function with
findgRNAs = FALSE
as described above.
inputFile1Path <- system.file("extdata", "rs362331C.fa", package="CRISPRseek")
inputFile2Path <- system.file("extdata", "rs362331T.fa", package="CRISPRseek")
res <- compare2Sequences(inputFile1Path,
inputFile2Path,
outputDir = outputDir,
overwrite = TRUE)
## search for gRNAs for input file1...
## search for gRNAs for input file2...
## [1] "Scoring ..."
## finish off-target search in sequence 2
## finish off-target search in sequence 1
## finish feature vector building
## finish score calculation
## [1] "Done!"
## name gRNAPlusPAM scoreForSeq1 scoreForSeq2
## 10 rs362331C_gR22r TGAAGTGCACACAGTGGATGNGG 100 17.2
## 11 rs362331C_gR21r GAAGTGCACACAGTGGATGANGG 100 26.8
## 12 rs362331C_gR15r CACACAGTGGATGAGGGAGCNGG 100 61.1
## 7 rs362331C_gR9r GTGGATGAGGGAGCAGGCGTNGG 100 98.6
## 2 rs362331T_gR38f CTACTGTGTGCACTTCATCCNGG 100 100
## 6 rs362331T_gR10r AGTAGATGAGGGAGCAGGCGNGG 100 100
## gRNAefficacy scoreDiff
## 10 0.0978516718170484 82.8
## 11 0.129816665166513 73.2
## 12 0.0169291602963948 38.9
## 7 extended sequence too short 1.4
## 2 extended sequence too short 0.0
## 6 0.161377523263354 0.0
Cytosine base editors (CBEs) and adenine base editors (ABEs) can
introduce specific DNA C-to-T or A-to-G alterations, respectively (Gaudelli et al. 2017;
Komor et al. 2016). The
offTargetAnalsys()
can design gRNAs optimized for base
editing by setting baseEditing = TRUE
. In this case, the
following parameters must also be specified (defaults are for the CBE
system developed in the Liu Lab): targetBase = "C"
,
editingWindow = 4:8
,
editingWindow.offtarget = 4:8
.
res <- offTargetAnalysis(inputFilePath = inputFilePath,
chromToSearch = "",
baseEditing = TRUE,
targetBase = "C",
editingWindow = 4:8,
editingWindow.offtargets = 4:8,
outputDir = outputDir,
overwrite = TRUE)
res
## DNAStringSet object of length 1:
## width seq names
## [1] 23 CCAGAGCAGGATCCACAAACTGG Hsap_GATA1_ex2_gR7r
In addition to CBE, the Liu Lab also developed the prime editor (PE) (Anzalone et al. 2019), which is more versatile and flexible with high efficacy and without the need to make a DSB or providing donor template. It can be used to make all possible 12 single base changes, 1-44 bp insertions, or 1-80 bp deletions. This editing system can be programmed to correct about 89 percent of human pathogenic variants.
To design gRNAs and pegRNAs for PE, set
primeEditing = TRUE
together with the following
parameters:
Type ?offTargetAnalysis
for detailed description of each
of these parameters.
inputFilePath <- DNAStringSet("CCAGTTTGTGGATCCTGCTCTGGTGTCCTCCACACCAGAATCAGGGATCGAAAACTCATCAGTCGATGCGAGTCATCTAAATTCCGATCAATTTCACACTTTAAACG")
res <- offTargetAnalysis(inputFilePath,
chromToSearch = "",
gRNAoutputName = "testPEgRNAs", # Required when inputFilePath is a DNAStringSet object.
primeEditing = TRUE,
targeted.seq.length.change = 0,
bp.after.target.end = 15,
PBS.length = 15,
RT.template.length = 8:30,
RT.template.pattern = "D$",
target.start = 20,
target.end = 20,
corrected.seq = "T",
findPairedgRNAOnly = TRUE,
paired.orientation = "PAMin",
outputDir = outputDir,
overwrite = TRUE)
res
## DNAStringSet object of length 4:
## width seq names
## [1] 23 CCAGTTTGTGGATCCTGCTCTGG NA_gR17f
## [2] 23 GAGTTTTCGATCCCTGATTCTGG NA_gR41r
## [3] 23 TTCGATCCCTGATTCTGGTGTGG NA_gR36r
## [4] 23 GATCCCTGATTCTGGTGTGGAGG NA_gR33r
For questions related to usage, please search/post your queries on the Bioconductor Support Site. If you wish to report a bug or request a new feature, kindly raise an issue on the CRISPRseek GitHub repository.
Yes, starting from version 1.44.1, CRISPRseek
supports the detection of off-targets with bulges by integrating
Cas-OFFinder (Bae, Park,
and Kim 2014). To learn more, type
?getOfftargetWithBulge
and ?offTargetAnalysis
for detailed documentation.
If you use CRISPRseek in your work, please cite it as follows:
## Please cite the paper below for the CRISPRseek package.
##
## Zhu LJ, Holmes BR, Aronin N, Brodsky MH (2014) CRISPRseek: A
## Bioconductor Package to Identify Target-Specific Guide RNAs for
## CRISPR-Cas9 Genome-Editing Systems. PLoS ONE 9(9): e108424.
## doi:10.1371/journal.pone.0108424
##
## Lihua Julie Zhu (2015). Overview of guide RNA design tools for
## CRISPR-Cas9 genome editing technology. Front. Biol., 10(4): 289-296
##
## To see these entries in BibTeX format, use 'print(<citation>,
## bibtex=TRUE)', 'toBibtex(.)', or set
## 'options(citation.bibtex.max=999)'.
Here is the output of sessionInfo()
on the system on
which this document was compiled running pandoc 3.2.1:
## R version 4.4.3 (2025-02-28)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.2 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] org.Hs.eg.db_3.20.0
## [2] TxDb.Hsapiens.UCSC.hg38.knownGene_3.20.0
## [3] BSgenome.Hsapiens.UCSC.hg38_1.4.5
## [4] BSgenome_1.75.1
## [5] rtracklayer_1.67.1
## [6] BiocIO_1.17.1
## [7] CRISPRseek_1.47.1
## [8] GenomicFeatures_1.59.1
## [9] AnnotationDbi_1.69.0
## [10] Biobase_2.67.0
## [11] GenomicRanges_1.59.1
## [12] Biostrings_2.75.4
## [13] GenomeInfoDb_1.43.4
## [14] XVector_0.47.2
## [15] IRanges_2.41.3
## [16] S4Vectors_0.45.4
## [17] BiocGenerics_0.53.6
## [18] generics_0.1.3
## [19] BiocFileCache_2.15.1
## [20] dbplyr_2.5.0
## [21] BiocStyle_2.35.0
##
## loaded via a namespace (and not attached):
## [1] DBI_1.2.3 bitops_1.0-9
## [3] rlang_1.1.5 magrittr_2.0.3
## [5] ade4_1.7-23 rio_1.2.3
## [7] matrixStats_1.5.0 compiler_4.4.3
## [9] RSQLite_2.3.9 png_0.1-8
## [11] vctrs_0.6.5 stringr_1.5.1
## [13] pkgconfig_2.0.3 crayon_1.5.3
## [15] fastmap_1.2.0 Rsamtools_2.23.1
## [17] rmarkdown_2.29 UCSC.utils_1.3.1
## [19] purrr_1.0.4 bit_4.6.0
## [21] xfun_0.51 cachem_1.1.0
## [23] seqinr_4.2-36 jsonlite_1.9.1
## [25] blob_1.2.4 rhdf5filters_1.19.2
## [27] keras_2.15.0 DelayedArray_0.33.6
## [29] Rhdf5lib_1.29.1 BiocParallel_1.41.2
## [31] parallel_4.4.3 tensorflow_2.16.0
## [33] R6_2.6.1 bslib_0.9.0
## [35] stringi_1.8.4 reticulate_1.41.0.1
## [37] jquerylib_0.1.4 Rcpp_1.0.14
## [39] SummarizedExperiment_1.37.0 knitr_1.49
## [41] base64enc_0.1-3 Matrix_1.7-2
## [43] tidyselect_1.2.1 abind_1.4-8
## [45] yaml_2.3.10 codetools_0.2-20
## [47] curl_6.2.1 lattice_0.22-6
## [49] tibble_3.2.1 withr_3.0.2
## [51] KEGGREST_1.47.0 evaluate_1.0.3
## [53] zip_2.3.2 pillar_1.10.1
## [55] BiocManager_1.30.25 filelock_1.0.3
## [57] MatrixGenerics_1.19.1 whisker_0.4.1
## [59] RCurl_1.98-1.16 gtools_3.9.5
## [61] glue_1.8.0 mltools_0.3.5
## [63] maketools_1.3.2 tools_4.4.3
## [65] sys_3.4.3 data.table_1.17.0
## [67] openxlsx_4.2.8 GenomicAlignments_1.43.0
## [69] buildtools_1.0.0 XML_3.99-0.18
## [71] rhdf5_2.51.2 grid_4.4.3
## [73] GenomeInfoDbData_1.2.13 restfulr_0.0.15
## [75] cli_3.6.4 tfruns_1.5.3
## [77] S4Arrays_1.7.3 dplyr_1.1.4
## [79] hash_2.2.6.3 zeallot_0.1.0
## [81] sass_0.4.9 digest_0.6.37
## [83] SparseArray_1.7.6 rjson_0.2.23
## [85] memoise_2.0.1 htmltools_0.5.8.1
## [87] lifecycle_1.0.4 httr_1.4.7
## [89] bit64_4.6.0-1 MASS_7.3-65