LRBase
and scTensor
scTensor 2.14.0
This vignette has been changed in BioC 3.14, when each data package (LRBase.XXX.eg.db) is deprecated and the way to provide LRBase data has changed to AnnotationHub-style.
LRBase
and scTensor
from BioC 3.14 (Nov. 2021)This section is for the users of previous LRBase.XXX.eg.db-type packages and scTensor. The specifications of the LRBase.XXX.eg.db and scTensor have changed significantly since BioC 3.14. Specifically, the distribution of all LRBase.XXX.eg.db-type packages will be abolished, and the policy has been switched to one where the data is placed on a cloud server called AnnotationHub, and users are allowed to retrieve the data only when they really need it. The following are the advantages of this AnnotationHub-style.
Due to the rapid development of single-cell RNA-Seq (scRNA-Seq) technologies, wide variety of cell types such as multiple organs of a healthy person, stem cell niche and cancer stem cell have been found. Such complex systems are composed of communication between cells (cell-cell interaction or CCI).
Many CCI studies are based on the ligand-receptor (L-R)-pair list of FANTOM5 project1 Jordan A. Ramilowski, A draft network of ligand-receptor-mediated multicellular signaling in human, Nature Communications, 2015 as the evidence of CCI (http://fantom.gsc.riken.jp/5/suppl/Ramilowski_et_al_2015/data/PairsLigRec.txt). The project proposed the L-R-candidate genes by following two basises.
The project also merged the data with previous L-R database such as IUPHAR/DLRP/HPMR and filter out the list without PMIDs. The recent L-R databases such as CellPhoneDB and SingleCellSignalR also manually curated L-R pairs, which are not listed in IUPHAR/DLRP/HPMR. In Bader Laboratory, many putative L-R databases are predicted by their standards. In our framework, we expanded such L-R databases for 134 organisms based on the ortholog relationships. For the details, check the summary of rikenbit/lrbase-workflow2 https://github.com/rikenbit/lrbase-workflow#summary, which is the Snakemake workflow to create LRBase data in each bi-annual update of Bioconductor.
LRBase
and scTensor
frameworkOur L-R databases (LRBase
) are provided a cloud server called AnnotationHub, and users are allowed to retrieve the data only when they really need it. Downloaded data is stored as a cache file on our local machines by the BiocFileCache mechanism. Then, the data is converted to LRBase object by LRBaseDbi. We also developed scTensor, which is a method to detect CCI and the CCI-related L-R pairs simultaneously. This document provides the way to use LRBaseDbi, LRBase objects, and scTensor (Figure 1).
To create the LRBase of 134 organisms, we introduced 36 approarches including known/putative L-R pairing. Please see the evidence code of lrbase-workflow3 https://github.com/rikenbit/lrbase-workflow.
AnnotationHub
First of all, we download the data of LRBase from AnnotationHub.
AnnotationHub::AnnotationHub
retrieve the metadata of all the data stored in cloud server.
library("AnnotationHub")
ah <- AnnotationHub()
mcols(ah)
## DataFrame with 71308 rows and 15 columns
## title dataprovider species taxonomyid
## <character> <character> <character> <integer>
## AH5012 Chromosome Band UCSC Homo sapiens 9606
## AH5013 STS Markers UCSC Homo sapiens 9606
## AH5014 FISH Clones UCSC Homo sapiens 9606
## AH5015 Recomb Rate UCSC Homo sapiens 9606
## AH5016 ENCODE Pilot UCSC Homo sapiens 9606
## ... ... ... ... ...
## AH116725 TENET_consensus_open.. NA Homo sapiens 9606
## AH116726 TENET_consensus_prom.. NA Homo sapiens 9606
## AH116727 ENCODE_dELS_regions ENCODE Homo sapiens 9606
## AH116728 ENCODE_pELS_regions ENCODE Homo sapiens 9606
## AH116729 ENCODE_PLS_regions ENCODE Homo sapiens 9606
## genome description coordinate_1_based
## <character> <character> <integer>
## AH5012 hg19 GRanges object from .. 1
## AH5013 hg19 GRanges object from .. 1
## AH5014 hg19 GRanges object from .. 1
## AH5015 hg19 GRanges object from .. 1
## AH5016 hg19 GRanges object from .. 1
## ... ... ... ...
## AH116725 hg38 A composite GRanges .. 1
## AH116726 hg38 A composite GRanges .. 1
## AH116727 hg38 A GRanges object con.. 1
## AH116728 hg38 A GRanges object con.. 1
## AH116729 hg38 A GRanges object con.. 1
## maintainer rdatadateadded preparerclass
## <character> <character> <character>
## AH5012 Marc Carlson <mcarls.. 2013-03-26 UCSCFullTrackImportP..
## AH5013 Marc Carlson <mcarls.. 2013-03-26 UCSCFullTrackImportP..
## AH5014 Marc Carlson <mcarls.. 2013-03-26 UCSCFullTrackImportP..
## AH5015 Marc Carlson <mcarls.. 2013-03-26 UCSCFullTrackImportP..
## AH5016 Marc Carlson <mcarls.. 2013-03-26 UCSCFullTrackImportP..
## ... ... ... ...
## AH116725 Rhie Lab at the Univ.. 2024-04-29 TENET.AnnotationHub
## AH116726 Rhie Lab at the Univ.. 2024-04-29 TENET.AnnotationHub
## AH116727 Rhie Lab at the Univ.. 2024-04-29 TENET.AnnotationHub
## AH116728 Rhie Lab at the Univ.. 2024-04-29 TENET.AnnotationHub
## AH116729 Rhie Lab at the Univ.. 2024-04-29 TENET.AnnotationHub
## tags rdataclass
## <AsIs> <character>
## AH5012 cytoBand,UCSC,track,... GRanges
## AH5013 stsMap,UCSC,track,... GRanges
## AH5014 fishClones,UCSC,track,... GRanges
## AH5015 recombRate,UCSC,track,... GRanges
## AH5016 encodeRegions,UCSC,track,... GRanges
## ... ... ...
## AH116725 ENCODE,GenomicSequence,Homo_sapiens,... GRanges
## AH116726 EpigenomeRoadMap,GenomicSequence,Homo_sapiens,... GRanges
## AH116727 DnaseSeq,ENCODE,GenomicSequence,... GRanges
## AH116728 DnaseSeq,ENCODE,GenomicSequence,... GRanges
## AH116729 DnaseSeq,ENCODE,GenomicSequence,... GRanges
## rdatapath sourceurl sourcetype
## <character> <character> <character>
## AH5012 goldenpath/hg19/data.. rtracklayer://hgdown.. UCSC track
## AH5013 goldenpath/hg19/data.. rtracklayer://hgdown.. UCSC track
## AH5014 goldenpath/hg19/data.. rtracklayer://hgdown.. UCSC track
## AH5015 goldenpath/hg19/data.. rtracklayer://hgdown.. UCSC track
## AH5016 goldenpath/hg19/data.. rtracklayer://hgdown.. UCSC track
## ... ... ... ...
## AH116725 TENET.AnnotationHub/.. https://github.com/r.. Multiple
## AH116726 TENET.AnnotationHub/.. https://github.com/r.. Multiple
## AH116727 TENET.AnnotationHub/.. https://screen.encod.. BED
## AH116728 TENET.AnnotationHub/.. https://screen.encod.. BED
## AH116729 TENET.AnnotationHub/.. https://screen.encod.. BED
Specifying some keywords in query()
, we can find LRBase data in AnnotationHub.
dbfile <- query(ah, c("LRBaseDb", "Homo sapiens", "v002"))[[1]]
## loading from cache
AnnotationHub also keeps old data as an archive, so please make sure you have the latest version (e.g. “v002” or higher) when you search for LRBaseDb.
Then, we can convert dbfile
into LRBase object by using LRBaseDbi
.
library("LRBaseDbi")
## LRBase.XXX.eg.db-type packages are deprecated since Bioconductor 3.14. Use AnnotationHub instead. For details, check the vignette of LRBaseDbi
LRBase.Hsa.eg.db <- LRBaseDbi::LRBaseDb(dbfile)
Some data access functions are available for LRBase objects.
Any data table are retrieved by 4 functions defined by
AnnotationDbi; columns
, keytypes
, keys
, and select
and commonly implemented by LRBaseDbi package. columns
returns the rows which we can retrieve in LRBase objects. keytypes
returns the rows which can be used as the optional parameter in keys
and select functions against LRBase objects. keys
function returns the value of keytype. select
function returns the rows in particular columns, which are having user-specified keys. This function returns the result as a dataframe. See the vignette of AnnotationDbi for more details.
columns(LRBase.Hsa.eg.db)
## [1] "GENEID_L" "GENEID_R" "SOURCEDB" "SOURCEID"
keytypes(LRBase.Hsa.eg.db)
## [1] "GENEID_L" "GENEID_R" "SOURCEDB" "SOURCEID"
key_HSA <- keys(LRBase.Hsa.eg.db, keytype="GENEID_L")
head(select(LRBase.Hsa.eg.db, keys=key_HSA[1:2],
columns=c("GENEID_L", "GENEID_R"), keytype="GENEID_L"))
## GENEID_L GENEID_R
## 1 1 6622
## 2 1 310
## 3 1 10321
## 4 1 10549
## 5 100 1803
Other additional functions like species
, nomenclature
, and listDatabases
are available. In each LRBase.XXX.eg.db-type package, species
function returns the common name and nomenclature
returns the scientific name. listDatabases
function returns the source of data. dbInfo
returns the information of the package. dbfile
returns the directory where sqlite file is stored. dbschema
returns the schema of the database. dbconn
returns the connection to the sqlite database.
lrNomenclature(LRBase.Hsa.eg.db)
## [1] "Homo sapiens"
species(LRBase.Hsa.eg.db)
## [1] "Human"
lrListDatabases(LRBase.Hsa.eg.db)
## SOURCEDB
## 1 BADERLAB
## 2 CELLPHONEDB
## 3 DLRP
## 4 FANTOM5
## 5 HPMR
## 6 IUPHAR
## 7 SINGLECELLSIGNALR
## 8 SOURCEDB
## 9 SWISSPROT_HPRD
## 10 SWISSPROT_STRING
## 11 TREMBL_HPRD
## 12 TREMBL_STRING
lrVersion(LRBase.Hsa.eg.db)
## NAME VALUE
## 1 LRVERSION 2018
dbInfo(LRBase.Hsa.eg.db)
## NAME
## 1 NAME
## 2 SOURCEDATE
## 3 SOURCENAME1
## 4 SOURCENAME2
## 5 SOURCENAME3
## 6 SOURCENAME4
## 7 SOURCENAME5
## 8 SOURCENAME6
## 9 SOURCENAME7
## 10 SOURCENAME8
## 11 SOURCENAME9
## 12 SOURCENAME10
## 13 SOURCENAME11
## 14 SOURCEURL1
## 15 SOURCEURL2
## 16 SOURCEURL3
## 17 SOURCEURL4
## 18 SOURCEURL5
## 19 SOURCEURL6
## 20 SOURCEURL7
## 21 SOURCEURL8
## 22 SOURCEURL9
## 23 SOURCEURL10
## 24 SOURCEURL11
## 25 DBSCHEMA
## 26 DBSCHEMAVERSION
## 27 ORGANISM
## 28 SPECIES
## 29 package
## 30 Db type
## 31 LRVERSION
## 32 TAXID
## VALUE
## 1 VALUE
## 2 2021-08-11
## 3 DLRP
## 4 IUPHAR
## 5 HPMR
## 6 CELLPHONEDB
## 7 SINGLECELLSIGNALR
## 8 FANTOM5
## 9 BADERLAB
## 10 SWISSPROT
## 11 HPRD
## 12 TREMBL
## 13 STRING
## 14 https://dip.doe-mbi.ucla.edu/dip/DLRP.cgi
## 15 https://www.guidetopharmacology.org/download.jsp
## 16
## 17 https://www.cellphonedb.org
## 18 http://www.bioconductor.org/packages/release/bioc/html/SingleCellSignalR.html
## 19 http://fantom.gsc.riken.jp
## 20 http://baderlab.org/CellCellInteractions
## 21 http://www.uniprot.org/uniprot/?query=reviewed:yes
## 22 http://hprd.org/download
## 23 http://www.uniprot.org/uniprot/?query=reviewed:no
## 24 https://string-db.org/cgi/download.pl
## 25 LRBase.Hsa.eg.db
## 26 1.0
## 27 Homo sapiens
## 28 Human
## 29 AnnotationDbi
## 30 LRBaseDb
## 31 2018
## 32 9606
dbfile(LRBase.Hsa.eg.db)
## AH97772
## "/home/biocbuild/.cache/R/AnnotationHub/2fc02a6abe6a9f_104518"
dbschema(LRBase.Hsa.eg.db)
## [1] "CREATE TABLE DATA (\n GENEID_L VARCHAR(10) NOT NULL, -- e.g., 19\n GENEID_R VARCHAR(10) NOT NULL, -- e.g., 3763409\n SOURCEID VARCHAR (10) NOT NULL, -- e.g., 27535533\n SOURCEDB VARCHAR (10) NOT NULL -- e.g., SWISSPROT_STRING\n)"
## [2] "CREATE TABLE METADATA (\n NAME NOT NULL,\n VALUE TEXT\n)"
## [3] "CREATE INDEX A ON DATA (GENEID_L)"
## [4] "CREATE INDEX B ON DATA (GENEID_R)"
## [5] "CREATE INDEX C ON DATA (SOURCEDB)"
## [6] "CREATE INDEX D ON DATA (SOURCEID)"
dbconn(LRBase.Hsa.eg.db)
## <SQLiteConnection>
## Path: /home/biocbuild/.cache/R/AnnotationHub/2fc02a6abe6a9f_104518
## Extensions: TRUE
Combined with dbGetQuery
function of RSQLite package,
more complicated queries also can be submitted.
suppressPackageStartupMessages(library("RSQLite"))
dbGetQuery(dbconn(LRBase.Hsa.eg.db),
"SELECT * FROM DATA WHERE GENEID_L = '9068' AND GENEID_R = '14' LIMIT 10")
## [1] GENEID_L GENEID_R SOURCEID SOURCEDB
## <0 rows> (or 0-length row.names)
scTensor
(CCI-tensor construction, decomposition, and HTML reporting)Combined with LRBase object and user’s gene expression matrix of scRNA-Seq, scTensor detects CCIs and generates HTML reports for exploratory data inspection. The algorithm of scTensor is as follows.
Firstly, scTensor calculates the celltype-level mean vectors, searches the corresponding pair of genes in the row names of the matrix, and extracted as two vectors.
Next, the cell type-level mean vectors of ligand expression and that of receptor expression are multiplied as outer product and converted to cell type \(\times\) cell type matrix. Here, the multiple matrices can be represented as a three-order “tensor” (Ligand-Cell * Receptor-Cell * L-R-Pair). scTensor decomposes the tensor into a small tensor (core tensor) and two factor matrices. Tensor decomposition is very similar to the matrix decomposition like PCA (principal component analysis). The core tensor is similar to the eigenvalue of PCA; this means that how much the pattern is outstanding. Likewise, three matrices are similar to the PC scores/loadings of PCA; These represent which ligand-cell/receptor-cell/L-R-pair are informative. When the matrices have negative values, interpreting which direction (+/-) is important and which is not, is a difficult and laboring task. That’s why, scTensor performs non-negative Tucker2 decomposition (NTD2), which is non-negative version of tensor decomposition (cf. nnTensor).
Finally, the result of NTD2 is summarized as an HTML report. Because most of the plots are visualized by plotly package, the precise information of the plot can be interactively confirmed by user’s on-site web browser. The two factor matrices can be interactively viewed and which cell types and which L-R-pairs are likely to be interacted each other. The mode-3 (LR-pair direction) sum of the core tensor is calculated and visualized as Ligand-Receptor Patterns. Detail of (Ligand-Cell, Receptor-Cell, L-R-pair) Patterns are also visualized.
SingleCellExperiment
objectHere, we use the scRNA-Seq dataset of male germline cells and somatic cells\(^{3}\)GSE86146 as demo data. For saving the package size, the number of genes is strictly reduced by the standard of highly variable genes with a threshold of the p-value are 1E-150 (cf. Identifying highly variable genes). That’s why we won’t argue about the scientific discussion of the data here.
We assume that user has a scRNA-Seq data matrix containing expression count data summarised at the level of the gene. First, we create a SingleCellExperiment object containing the data. The rows of the object correspond to features, and the columns correspond to cells. The gene identifier is limited as NCBI Gene ID for now.
To improve the interpretability of the following HTML report, we highly recommend that user specifies the two-dimensional data of input data (e.g. PCA, t-SNE, or UMAP). Such information is easily specified by reducedDims
function of SingleCellExperiment package and is saved to reducedDims slot of SingleCellExperiment
object (Figure ??).
suppressPackageStartupMessages(library("scTensor"))
suppressPackageStartupMessages(library("SingleCellExperiment"))
data(GermMale)
data(labelGermMale)
data(tsneGermMale)
sce <- SingleCellExperiment(assays=list(counts = GermMale))
reducedDims(sce) <- SimpleList(TSNE=tsneGermMale$Y)
plot(reducedDims(sce)[[1]], col=labelGermMale, pch=16, cex=2,
xlab="Dim1", ylab="Dim2", main="Germline, Male, GSE86146")
legend("topleft", legend=c(paste0("FGC_", 1:3), paste0("Soma_", 1:4)),
col=c("#9E0142", "#D53E4F", "#F46D43", "#ABDDA4", "#66C2A5", "#3288BD", "#5E4FA2"),
pch=16)
cellCellSetting
To perform the tensor decomposition and HTML report, user is supposed to specify
to SingleCellExperiment
object.
cellCellSetting(sce, LRBase.Hsa.eg.db, names(labelGermMale))
The corresponding information is registered to the metadata slot of SingleCellExperiment
object by cellCellSetting
function.
cellCellDecomp
After cellCellSetting
, we can perform tensor decomposition by cellCellDecomp
. Here the parameter ranks
is specified as dimension of core tensor. For example, c(2, 3) means The data tensor is decomposed to 2 ligand-patterns and 3 receptor-patterns.
set.seed(1234)
cellCellDecomp(sce, ranks=c(2,3))
## Input data matrix may contains 7 gene symbols because the name contains some alphabets.
## scTensor uses only NCBI Gene IDs for now.
## Here, the gene symbols are removed and remaining 235 NCBI Gene IDs are used for scTensor next step.
## 7 * 7 * 4 Tensor is created
Although user has to specify the rank to perform cellCellDecomp, we implemented a simple rank estimation function based on the eigenvalues distribution of PCA in the matricised tensor in each mode in cellCellRank
. rks$selected
is also specified as rank parameter of cellCellDecomp
.
(rks <- cellCellRanks(sce))
## Each rank, multiple NMF runs are performed
##
|
| | 0%
|
|========== | 14%
|
|==================== | 29%
|
|============================== | 43%
|
|======================================== | 57%
|
|================================================== | 71%
|
|============================================================ | 86%
|
|======================================================================| 100%
## Each rank estimation method
##
|
| | 0%
|
|========== | 14%
|
|==================== | 29%
|
|============================== | 43%
|
|======================================== | 57%
|
|================================================== | 71%
|
|============================================================ | 86%
|
|======================================================================| 100%
## Each rank, multiple NMF runs are performed
##
|
| | 0%
|
|========== | 14%
|
|==================== | 29%
|
|============================== | 43%
|
|======================================== | 57%
|
|================================================== | 71%
|
|============================================================ | 86%
|
|======================================================================| 100%
## Each rank estimation method
##
|
| | 0%
|
|========== | 14%
|
|==================== | 29%
|
|============================== | 43%
|
|======================================== | 57%
|
|================================================== | 71%
|
|============================================================ | 86%
|
|======================================================================| 100%
## $RSS
## $RSS$rss1
## [1] 124587604751 25929822763 2085959976 438461900 312129322
## [6] 277726907 183431820
##
## $RSS$rss2
## [1] 122493907374 6776993503 2568919182 3253818525 1438749444
## [6] 414514697 1111518908
##
##
## $selected
## [1] 3 2
rks$selected
## [1] 3 2
cellCellReport
If cellCellDecomp
is properly finished, we can perform cellCellReport
function to output the HTML report like below. Please type example(cellCellReport)
and the report will be generated in the temporary directory (it costs 5 to 10 minutes). After cellCellReport
, multiple R markdown files, compiled HTML files, figures, and R binary file containing the result of analysis are saved to out.dir
(Figure 2). For more details, open the index.html
by your web browser. Combined with cloud storage service such as Amazon Simple Storage Service (S3), it can be a simple web application and multiple people like collaborators can confirm the same report simultaneously.
## R version 4.4.0 beta (2024-04-15 r86425)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 22.04.4 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.19-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_US.UTF-8 LC_COLLATE=en_US.UTF-8
## [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] SingleCellExperiment_1.26.0
## [2] SummarizedExperiment_1.34.0
## [3] MatrixGenerics_1.16.0
## [4] matrixStats_1.3.0
## [5] RSQLite_2.3.6
## [6] LRBaseDbi_2.14.0
## [7] AnnotationHub_3.12.0
## [8] BiocFileCache_2.12.0
## [9] dbplyr_2.5.0
## [10] scTGIF_1.18.0
## [11] Homo.sapiens_1.3.1
## [12] TxDb.Hsapiens.UCSC.hg19.knownGene_3.2.2
## [13] org.Hs.eg.db_3.19.1
## [14] GO.db_3.19.1
## [15] OrganismDbi_1.46.0
## [16] GenomicFeatures_1.56.0
## [17] GenomicRanges_1.56.0
## [18] GenomeInfoDb_1.40.0
## [19] AnnotationDbi_1.66.0
## [20] IRanges_2.38.0
## [21] S4Vectors_0.42.0
## [22] Biobase_2.64.0
## [23] BiocGenerics_0.50.0
## [24] scTensor_2.14.0
## [25] BiocStyle_2.32.0
##
## loaded via a namespace (and not attached):
## [1] fs_1.6.4 bitops_1.0-7 enrichplot_1.24.0
## [4] HDO.db_0.99.1 httr_1.4.7 webshot_0.5.5
## [7] RColorBrewer_1.1-3 Rgraphviz_2.48.0 tools_4.4.0
## [10] backports_1.4.1 utf8_1.2.4 R6_2.5.1
## [13] lazyeval_0.2.2 withr_3.0.0 prettyunits_1.2.0
## [16] graphite_1.50.0 gridExtra_2.3 schex_1.18.0
## [19] fdrtool_1.2.17 cli_3.6.2 TSP_1.2-4
## [22] scatterpie_0.2.2 entropy_1.3.1 sass_0.4.9
## [25] genefilter_1.86.0 meshr_2.10.0 Rsamtools_2.20.0
## [28] yulab.utils_0.1.4 gson_0.1.0 txdbmaker_1.0.0
## [31] DOSE_3.30.0 MeSHDbi_1.40.0 AnnotationForge_1.46.0
## [34] nnTensor_1.2.0 plotrix_3.8-4 maps_3.4.2
## [37] visNetwork_2.1.2 generics_0.1.3 gridGraphics_0.5-1
## [40] GOstats_2.70.0 BiocIO_1.14.0 dplyr_1.1.4
## [43] dendextend_1.17.1 Matrix_1.7-0 fansi_1.0.6
## [46] abind_1.4-5 lifecycle_1.0.4 yaml_2.3.8
## [49] qvalue_2.36.0 SparseArray_1.4.0 grid_4.4.0
## [52] blob_1.2.4 misc3d_0.9-1 crayon_1.5.2
## [55] lattice_0.22-6 msigdbr_7.5.1 cowplot_1.1.3
## [58] annotate_1.82.0 KEGGREST_1.44.0 magick_2.8.3
## [61] pillar_1.9.0 knitr_1.46 fgsea_1.30.0
## [64] tcltk_4.4.0 rjson_0.2.21 codetools_0.2-20
## [67] fastmatch_1.1-4 glue_1.7.0 outliers_0.15
## [70] ggfun_0.1.4 data.table_1.15.4 vctrs_0.6.5
## [73] png_0.1-8 treeio_1.28.0 spam_2.10-0
## [76] rTensor_1.4.8 gtable_0.3.5 assertthat_0.2.1
## [79] cachem_1.0.8 xfun_0.43 mime_0.12
## [82] S4Arrays_1.4.0 tidygraph_1.3.1 survival_3.6-4
## [85] seriation_1.5.5 iterators_1.0.14 tinytex_0.50
## [88] fields_15.2 nlme_3.1-164 Category_2.70.0
## [91] ggtree_3.12.0 bit64_4.0.5 progress_1.2.3
## [94] filelock_1.0.3 bslib_0.7.0 colorspace_2.1-0
## [97] DBI_1.2.2 tidyselect_1.2.1 bit_4.0.5
## [100] compiler_4.4.0 curl_5.2.1 httr2_1.0.1
## [103] graph_1.82.0 xml2_1.3.6 DelayedArray_0.30.0
## [106] plotly_4.10.4 bookdown_0.39 shadowtext_0.1.3
## [109] rtracklayer_1.64.0 checkmate_2.3.1 scales_1.3.0
## [112] hexbin_1.28.3 RBGL_1.80.0 plot3D_1.4.1
## [115] rappdirs_0.3.3 stringr_1.5.1 digest_0.6.35
## [118] rmarkdown_2.26 ca_0.71.1 XVector_0.44.0
## [121] htmltools_0.5.8.1 pkgconfig_2.0.3 highr_0.10
## [124] fastmap_1.1.1 rlang_1.1.3 htmlwidgets_1.6.4
## [127] UCSC.utils_1.0.0 farver_2.1.1 jquerylib_0.1.4
## [130] jsonlite_1.8.8 BiocParallel_1.38.0 GOSemSim_2.30.0
## [133] RCurl_1.98-1.14 magrittr_2.0.3 GenomeInfoDbData_1.2.12
## [136] ggplotify_0.1.2 dotCall64_1.1-1 patchwork_1.2.0
## [139] munsell_0.5.1 Rcpp_1.0.12 babelgene_22.9
## [142] ape_5.8 viridis_0.6.5 stringi_1.8.3
## [145] tagcloud_0.6 ggraph_2.2.1 zlibbioc_1.50.0
## [148] MASS_7.3-60.2 plyr_1.8.9 parallel_4.4.0
## [151] ggrepel_0.9.5 Biostrings_2.72.0 graphlayouts_1.1.1
## [154] splines_4.4.0 hms_1.1.3 igraph_2.0.3
## [157] biomaRt_2.60.0 reshape2_1.4.4 BiocVersion_3.19.1
## [160] XML_3.99-0.16.1 evaluate_0.23 BiocManager_1.30.22
## [163] foreach_1.5.2 tweenr_2.0.3 tidyr_1.3.1
## [166] purrr_1.0.2 polyclip_1.10-6 heatmaply_1.5.0
## [169] ggplot2_3.5.1 ReactomePA_1.48.0 ggforce_0.4.2
## [172] xtable_1.8-4 restfulr_0.0.15 reactome.db_1.88.0
## [175] tidytree_0.4.6 viridisLite_0.4.2 tibble_3.2.1
## [178] aplot_0.2.2 ccTensor_1.0.2 memoise_2.0.1
## [181] registry_0.5-1 GenomicAlignments_1.40.0 cluster_2.1.6
## [184] concaveman_1.1.0 GSEABase_1.66.0