vignettes/sSNAPPY.Rmd
sSNAPPY.RmdThis vignette demonstrates how to use the package
sSNAPPY to compute directional single sample pathway
perturbation scores by incorporating pathway topologies, utilize sample
permutation to test the significance of individual scores and compare
average pathway activities across treatments.
The package also provides a function to visualise overlap between pathway genes contained in perturbed biological pathways as network plots.
The package sSNAPPY can be installed using the package
BiocManager, along with other packages required for this
vignette (tidyverse, magrittr,
ggplot2, cowplot, and DT).
if (!"BiocManager" %in% rownames(installed.packages()))
install.packages("BiocManager")
pkg <- c("tidyverse", "magrittr", "ggplot2", "cowplot", "DT")
BiocManager::install(pkg)
BiocManager::install("sSNAPPY")The example dataset used for this tutorial can loaded with
data() as shown below. It’s a subset of data retrieved from
Singhal
H et al. 2016, where ER-positive primary breast cancer tumour tissues
collected from 12 patients were split into fragments of equal sizes for
different treatments. For the purpose of this tutorial, we only took the
RNA-seq data from samples that were treated with vehicle, E2 OR E2 +
R5020 for 48 hrs. They were from 5 different patients, giving rise to 15
samples in total. More detailed description of the dataset can be
assessed through the help page (?logCPM_example and
?metadata_example).
data(logCPM_example)
data(metadata_example)
# check if samples included in the logCPM matrix and metadata dataframe are identical
setequal(colnames(logCPM_example), metadata_example$sample)## [1] TRUE
# View sample metadata
datatable(metadata_example, filter = "top")It is expected that the logCPM matrix will be filtered to remove undetectable genes and normalised to correct for library sizes or other systematic artefacts, such as gene length or GC contents, prior to applying this method. Filtration and normalisation has been performed on the example dataset.
Before single-sample logFCs (ssLogFCs) can be computed, row names of
the logCPM matrix need to be converted to entrez ID. This
is because all the pathway topology information retrieved will be in
entrez ID. The conversion can be achieved through
bioconductor packages AnnotationHub and
ensembldb.
head(logCPM_example)## Vehicle_N2_48 E2+R5020_N2_48 R5020_N2_48 Vehicle_N3_48 E2+R5020_N3_48
## 7105 4.081669 5.109147 4.695391 4.673548 5.217774
## 57147 5.211258 3.844185 4.343678 3.601125 3.427482
## 55732 2.111865 1.366069 2.471200 1.975589 1.713835
## 2268 6.930587 4.341791 4.388130 5.309135 5.066294
## 90529 1.725384 2.997631 3.035200 3.161472 3.320000
## 22875 4.099909 3.635736 3.845961 4.701901 4.886946
## R5020_N3_48 Vehicle_P4_48 E2+R5020_P4_48 R5020_P4_48 Vehicle_P5_48
## 7105 4.296947 4.667555 4.741936 4.621193 5.484910
## 57147 4.277160 5.128271 4.486209 4.624768 3.549639
## 55732 2.222948 2.437893 3.147235 3.327889 2.368816
## 2268 3.775364 2.912324 3.911222 4.258518 3.029545
## 90529 2.719585 4.578478 3.624418 3.660976 3.477681
## 22875 4.624176 5.453882 4.616053 5.200638 4.080996
## E2+R5020_P5_48 R5020_P5_48 Vehicle_P6_48 E2+R5020_P6_48 R5020_P6_48
## 7105 5.678579 5.118092 5.636142 5.434636 5.753488
## 57147 4.113434 3.802290 3.815933 3.965867 3.914067
## 55732 2.438389 2.594323 1.028753 1.493969 2.365653
## 2268 1.530079 3.044646 2.904204 3.796048 3.336565
## 90529 3.307986 3.073607 3.577813 3.450174 3.613297
## 22875 3.696430 3.267462 4.436550 4.363616 4.402518
To compute the ssLogFCs, samples must be in matching pairs. In the
example, treated samples are matched to the corresponding control sample
that were derived from the same patients. So the factor
parameter of the weight_ss_fc() functions needs to be set
to be “patient”. The function also requires the control treatment level
to be specified, which was “Vehicle” in this case.
weight_ss_fc() requires both the logCPM matrix and
sample metadata as input. The column names of the logCPM matrix should
be sample name, matching to a column called sample in the
metadata. The metadata must also contain a treatment column, and a
column corresponding to the factor parameter (i.e.. patient
in this case).
#compute weighted single sample logFCs
weightedFC <- weight_ss_fc(logCPM_example, metadata = metadata_example,
factor = "patient", control = "Vehicle")The weight_ss_fc() function firstly computes raw
ssLogFCs for each gene by subtracting logCPM values of control sample
from the logCPM values of treated samples for each patient.
It has been demonstrated previously that in RNA-seq data, lowly
expressed genes turn to have larger variance, which is also demonstrated
by the plots below. To reduce the impact of this artefact,
weight_ss_fc also weight each ssLogFCs by estimating the
relationship between the variance in ssLogFCs and mean logCPM, and
defining the gene-wise weight to be the inverse of the predicted
variance of that gene’s mean logCPM value.
perSample_FC <- lapply(levels(metadata_example$patient), function(x){
temp <- logCPM_example[seq_len(1000),str_detect(colnames(logCPM_example), x)]
ratio <- temp[, str_detect(colnames(temp), "Vehicle", negate = TRUE)] - temp[, str_detect(colnames(temp), "Vehicle")]
}) %>%
do.call(cbind,.)
aveCPM <- logCPM_example[seq_len(1000),] %>%
rowMeans() %>%
enframe(name = "gene_id",
value = "aveCPM")
p1 <- perSample_FC %>%
as.data.frame() %>%
rownames_to_column("gene_id") %>%
pivot_longer(cols = -"gene_id",
names_to = "name",
values_to = "logFC") %>%
left_join(aveCPM) %>%
ggplot(aes(aveCPM, logFC)) +
geom_point() +
labs(y = "sslogFC",
x = "Average logCPM") +
theme(
panel.background = element_blank()
)
p2 <- data.frame(
gene_id = rownames(perSample_FC),
variance = perSample_FC %>%
apply(1,var)) %>%
left_join(aveCPM) %>%
ggplot(aes(aveCPM, variance)) +
geom_point() +
geom_smooth(method = "loess") +
labs(y = "Variance in ssLogFCs",
x = "Average logCPM") +
theme(
panel.background = element_blank()
)
plot_grid(p1, p2)
The output of the weight_ss_fc() function is a list with
two element, where one is the weighted ssLogFCs matrix and the other is
a vector of gene-wise weights.
sSNAPPY adopts the pathway perturbation scoring algorithm proposed in SPIA, which makes use of gene-set topologies and gene-gene interaction to propagate pathway genes’ logFCs down the topologies to compute pathway perturbation scores, where signs of scores reflect the potential directions of changes.
Therefore, pathway topology information need to be firstly retrieved
from your chosen database and converted to weight adjacency matrices,
the format required to apply the scoring algorithm. This step is
achieved through a chain of functions that are part of the grapghite
and has been nested into one simple function in this package:
retrieve_topology(). Databases that are currently supported
are:
library(graphite)
graphite::pathwayDatabases() %>%
dplyr::filter(species == "hsapiens") %>%
pander::pander()| species | database |
|---|---|
| hsapiens | kegg |
| hsapiens | panther |
| hsapiens | pathbank |
| hsapiens | pharmgkb |
| hsapiens | reactome |
| hsapiens | smpdb |
| hsapiens | wikipathways |
The retrieved topology information will be a list where each element is a pathway. We recommend you to save it as a file so this step only needs to be performed once for each database.
This vignette chose KEGG pathways as an example. To reduce the computation time in the following steps, only half of the randomly sampled KEGG pathways were kept.
gsTopology <- retrieve_topology(database = "kegg")
head(names(gsTopology))## [1] "Glycolysis / Gluconeogenesis"
## [2] "Citrate cycle (TCA cycle)"
## [3] "Pentose phosphate pathway"
## [4] "Pentose and glucuronate interconversions"
## [5] "Fructose and mannose metabolism"
## [6] "Galactose metabolism"
If only selected pathways are of interest, it’s possible to only retrieve the topologies of those pathways by specifying the pathway names.
gsTopology_sub <- retrieve_topology(
database = "kegg",
pathwayName = c(
"Glycolysis / Gluconeogenesis",
"Citrate cycle (TCA cycle)",
"Pentose phosphate pathway"
))
names(gsTopology_sub)## [1] "Glycolysis / Gluconeogenesis" "Citrate cycle (TCA cycle)"
## [3] "Pentose phosphate pathway"
Once the expression matrix, sample metadata and pathway topologies are all ready, gene-wise single-sample perturbation scores can be computed:
genePertScore <- raw_gene_pert(weightedFC$logFC, gsTopology)and summed to derive pathway perturbation scores for each pathway in each treated samples.
pathwayPertScore <- pathway_pert(genePertScore)
head(pathwayPertScore)## sample tA gs_name
## 1 E2+R5020_N2_48 0.006751152 EGFR tyrosine kinase inhibitor resistance
## 2 R5020_N2_48 0.006744556 EGFR tyrosine kinase inhibitor resistance
## 3 E2+R5020_N3_48 0.002683974 EGFR tyrosine kinase inhibitor resistance
## 4 R5020_N3_48 -0.003974106 EGFR tyrosine kinase inhibitor resistance
## 5 E2+R5020_P4_48 0.004130870 EGFR tyrosine kinase inhibitor resistance
## 6 R5020_P4_48 0.015279018 EGFR tyrosine kinase inhibitor resistance
To derive the empirical p-values for each single sample PS or normalize the raw scores for comparing overall treatment effects, null distributions of scores for each pathway is generated through a sample-label permutation approach.
For each round of permutation, sample labels are randomly shuffled to
derive the permuted ssLogFCs, which are then used to score pathway
perturbation. We recommend to perform a minimum of 1000 rounds of
permutation, which means at least 8 samples are required. The
generate_permuted_scores() function does not require sample
metadata but the number of treatments in the study design, including the
control treatment, need to be specified. In this example data, the
number of treatment was 3.
Output of the generate_permuted_scores() function is a
list where each element is a vector of permuted perturbation scores for
a pathway.
The permutation step relies on the parallel computing feature
provided by BiocParallel.
User can choose to customize the parallel back-end or stick with the
default one returned by BiocParallel::bpparam(). Depending
on the size of the data, this step can take some time to complete. If
the sample size is large, we recommend users to consider performing this
step on a HPC.
permutedScore <- generate_permuted_scores(
expreMatrix = logCPM_example,
numOfTreat = 3, NB = 1000,
gsTopology = gsTopology,
weight = weightedFC$weight
)After the empirical null distributions are generated, the median and
mad will be calculated for each pathway to convert the test
single-sample perturbation scores derived from the
compute_perturbation_score() to robust z-scores: \[ (Score - Median)/MAD\] Two-sided p-values
associated with each robust z-scores are also computed and will be
corrected for multiple-testing using a user-define approach. The default
is fdr.
The normalise_by_permu function requires the test
perturbation score and permuted perturbation scores as input.
normalisedScores <- normalise_by_permu(permutedScore, ssPertScore)Since the permutation step takes a long time to run and the output is
too large to be included as part of the package, results of the
normalise_by_permu step has been pre-computed and can be
loaded with:
load(system.file("extdata", "normalisedScores.rda", package = "sSNAPPY"))The pathways that were significant perturbed within individual samples are:
normalisedScores %>%
dplyr::filter(adjPvalue < 0.05) %>%
left_join(metadata_example) %>%
mutate_at(vars(c("sample", "gs_name")), as.factor) %>%
mutate_if(is.numeric, sprintf, fmt = '%#.4f') %>%
mutate(Direction = ifelse(robustZ < 0, "Inhibited", "Activation")) %>%
dplyr::select(
sample, patient, Treatment = treatment,
`Perturbation Score` = robustZ, Direction,
`Gene-set name` = gs_name,
`P-value` = pvalue,
FDR = adjPvalue
) %>%
datatable(
filter = "top",
options = list(
columnDefs = list(list(targets = "Direction", visible = FALSE))
)
) %>%
formatStyle(
'Perturbation Score', 'Direction',
backgroundColor = styleEqual(c("Inhibited", "Activation"), c('lightblue', 'indianred'))
)We can use the plot_gs_network function to visualise
significantly perturbed biological pathways as networks, where edges
between gene-sets reflect how much overlap those two gene-sets share.
The function can take normalise_by_permu’s output, or a
subset of it as its direct input.
Nodes in the network plots could be coloured by the predicted direction of perturbation (i.e.. sign of robust z-score):
pl <- normalisedScores %>%
dplyr::filter(adjPvalue < 0.05) %>%
split(f = .$sample) %>%
lapply(
plot_gs_network,
# layout = "dh",
gsTopology = gsTopology,
colorBy = "robustZ"
) %>%
lapply(function(x){
x + theme(
panel.grid = element_blank(),
panel.background = element_blank()
) })
plot_grid(
plotlist = pl,
labels = names(pl),
label_size = 8,
nrow = 1)
Pathways significantly perturbed in individual samples, where gene-sets were colored by pathways’ direction of genes.
Or p-values:
pl <- normalisedScores %>%
dplyr::filter(adjPvalue < 0.05) %>%
split(f = .$sample) %>%
lapply(
plot_gs_network,
gsTopology = gsTopology,
colorBy = "pvalue",
color_lg_title = "P-value"
) %>%
lapply(function(x){
x + theme(
panel.grid = element_blank(),
panel.background = element_blank()
) })
plot_grid(
plotlist = pl,
labels = names(pl),
label_size = 8,
nrow = 1)
Pathways significantly perturbed in individual samples, where gene-sets were colored by pathways’ p-values
The function allows you to customize the layout, colour, edge
transparency and other aesthetics of the graph. More information can be
found in the help page (?plot_gs_network). Output of the
graph is a ggplot object and the theme of it can be changed
just as other ggplot figures.
In addtion to testing pathway perturbation at single-sample level, normalised perturbation scores can also be used to model mean treatment effects within a group, which in this case is the treatments. An advantage of this method is that it has a high level of flexibility that allows us to incorporate confounding factors cofactors or covariate to offset their effects.
For example, in the example data-set, samples were collected from patients with different progesterone receptor (PR) status and we can include it as a cofactor to offset its confounding effect.
fit <- normalisedScores %>%
left_join(metadata_example) %>%
split(f = .$gs_name) %>%
#.["Estrogen signaling pathway"] %>%
lapply(function(x)lm(robustZ ~ 0 + treatment + PR, data = x)) %>%
lapply(summary)
treat_sig <- lapply(
names(fit),
function(x){
fit[[x]]$coefficients %>%
as.data.frame() %>%
.[seq_len(2),] %>%
dplyr::select(Estimate, pvalue = `Pr(>|t|)` ) %>%
rownames_to_column("Treatment") %>%
mutate(
gs_name = x,
FDR = p.adjust(pvalue, "fdr"),
Treatment = str_remove_all(Treatment, "treatment")
)
}) %>%
bind_rows() Results from the linear modelling revealed pathways that were on average perturbed due to each treatment:
treat_sig %>%
dplyr::filter(FDR < 0.05) %>%
mutate_at(vars(c("Treatment", "gs_name")), as.factor) %>%
mutate_if(is.numeric, sprintf, fmt = '%#.4f') %>%
mutate(Direction = ifelse(Estimate < 0, "Inhibited", "Activation")) %>%
dplyr::select(
Treatment, `Perturbation Score` = Estimate, Direction,
`Gene-set name` = gs_name,
`P-value` = pvalue,
FDR
) %>%
datatable(
filter = "top",
options = list(
columnDefs = list(list(targets = "Direction", visible = FALSE))
)
) %>%
formatStyle(
'Perturbation Score', 'Direction',
backgroundColor = styleEqual(c("Inhibited", "Activation"), c('lightblue', 'indianred'))
)It is not surprising to see that the estrogen signalling pathway was significantly activated under both E2 and E2+R5020 treatments.
If there’s a specific pathway that we would like to dig deeper into
explore which pathway genes played a key role in the perturbation, for
example, activation of the “Ovarian steroidogenesis”, we can plot the
gene-level perturbation scores of all genes contained in that pathway as
a heatmap and annotate each column (ie. each sample) by the direction of
pathway perturbation in that sample or any other attributes using the
plot_gene_contribution() function.
We can see from the heatmap below that the activation of pathway was driven by a few dominating genes.
plot_gene_contribution(
genePertScore = genePertScore, gsToPlot = "Ovarian steroidogenesis", metadata = metadata_example,
annotation_attribute = c("pathwayPertScore", "treatment"), pathwayPertScore = pathwayPertScore,
breaks = seq(-0.0001, 0.0001, length.out = 100),
show_rownames = FALSE
)
Gene-level perturbation scores of all genes in the Ovarian steroidogenesis gene-set. The activation of pathway was driven by a few dominating genes.
By default, genes’ entrez ID are used and plotted as rownames, which
may not be very informative but the rownames could be overwritten with
the mapRownameTo parameter. Mapping between different gene
identifiers could be achieved through the mapIDs() function
from the Bioconductor package AnnotationDbi.
But to reduce the compiling time of this vignette, mapping between
entrez ID and gene names has been provided as a data.frame
called entrez2name.
Note that if the mapping information was provided and the mapping was successful for some genes but not the others, only genes that have been mapped to another identifier will be plotted.
load(system.file("extdata", "entrez2name.rda", package = "sSNAPPY"))
plot_gene_contribution(
genePertScore = genePertScore, gsToPlot = "Ovarian steroidogenesis", metadata = metadata_example,
annotation_attribute = c("pathwayPertScore", "treatment"), pathwayPertScore = pathwayPertScore,
breaks = seq(-0.0001, 0.0001, length.out = 100),
fontsize_row = 6,
mapEntrezID = entrez2name
)
Gene-level perturbation scores of all genes in the PI3K-Akt signaling pathway gene-set. A list of genes predominantly promoting the activation of the pathways could be identified.
Results of lm can also be visualised using the
plot_gs_network() function. We just need to change the name
of the Estimate column to robustZ to colour
the networks by the predicted directionality.
treat_sig %>%
dplyr::filter(FDR < 0.05, Treatment == "E2+R5020") %>%
dplyr::rename(robustZ = Estimate) %>%
plot_gs_network(
gsTopology = gsTopology,
colorBy = "robustZ"
) +
theme(
panel.grid = element_blank(),
panel.background = element_blank()
)
Pathways significantly perturbed by the E2+R5020 combintation treatment, where colors of nodes reflect pathways’ directions of changes.
By default, plot_gs_network() function does not include
nodes that are not connected to any other nodes, which could be turnt of
by setting the plotIsolated to TURE.
When a large number of pathways were perturbed, it is hard to answer
the question “What key biological processes were perturbed?”. To solve
that, we can use the plot_community() function to apply a
community detection algorithm to the network we had above and annotate
each community by the biological process that most pathways assigned to
that community belong to.
treat_sig %>%
dplyr::filter(FDR < 0.05, Treatment == "E2+R5020") %>%
dplyr::rename(robustZ = Estimate) %>%
plot_community(
gsTopology = gsTopology,
colorBy = "community",
color_lg_title = "Community"
) +
theme(panel.background = element_blank())
Pathways significantly perturbed by the E2+R5020 combintation treatment, annotated by the main biological processes each pathways belong to.
If we want to not only know if two pathways are connected but also
the genes connecting those pathways, we can use the
plot_gs2gene() function instead:
treat_sig %>%
dplyr::filter(FDR < 0.05, Treatment == "E2+R5020") %>%
dplyr::rename(robustZ = Estimate) %>%
plot_gs2gene(
gsTopology = gsTopology,
colorGS_By = "robustZ",
label_Gene = FALSE,
GeneNode_size = 1,
edgeAlpha = 0.1
) +
theme(panel.background = element_blank())
Pathways significantly perturbed by the E2+R5020 combintation treatment and all genes contained in those pathways.
However, since there are a large number of genes in each pathway, the plot above was quite message and it was unrealistic to plot all gene names. So it is recommend to filter pathway genes by their perturbation scores or logFCs.
For example, we can rank genes by the absolute values of their mean single-sample logFCs and only focus on genes that were ranked in the top 200.
meanFC <- weightedFC$logFC %>%
.[, str_detect(colnames(.), "E2", negate = TRUE)] %>%
apply(1, mean )
top200_gene <- meanFC %>%
abs() %>%
sort(decreasing = TRUE, ) %>%
.[1:200] %>%
names()
top200_FC <- meanFC %>%
.[names(.) %in% top200_gene]When we genes’ logFCs as a named vector, only pathway genes with
logFCs provided will be plotted and gene nodes will be colored by genes’
directions of changes. The names of logFC vector must be genes’ entrez
ID in the format of “ENTREZID:XXXX”, as pathway topology matrices
retrieved through retrieve_topology() always use entrez ID
as identifiers.
However, it is not going to be informative to label genes with their
entrez ID in the plot. So the same as the
plot_gene_contribution() function, we can provide a
data.frame to match genes’ entrez ID to our chosen
identifiers through the mapEntrezID parameter.
treat_sig %>%
dplyr::filter(FDR < 0.05, Treatment == "E2+R5020") %>%
dplyr::rename(robustZ = Estimate) %>%
plot_gs2gene(
gsTopology = gsTopology,
colorGS_By = "robustZ",
mapEntrezID = entrez2name,
geneFC = top200_FC,
edgeAlpha = 0.3,
GsName_size = 4
) +
theme(panel.background = element_blank())
Pathways significantly perturbed by the E2+R5020 combintation treatment and pathway genes with top 200 magnitudes of changes among al E2+R5020-treated samples, where both pathways and genes were colored by their directions of changes.
## R version 4.2.2 (2022-10-31)
## Platform: aarch64-apple-darwin20 (64-bit)
## Running under: macOS Monterey 12.2.1
##
## Matrix products: default
## BLAS: /Library/Frameworks/R.framework/Versions/4.2-arm64/Resources/lib/libRblas.0.dylib
## LAPACK: /Library/Frameworks/R.framework/Versions/4.2-arm64/Resources/lib/libRlapack.dylib
##
## locale:
## [1] en_AU.UTF-8/en_AU.UTF-8/en_AU.UTF-8/C/en_AU.UTF-8/en_AU.UTF-8
##
## attached base packages:
## [1] stats graphics grDevices utils datasets methods base
##
## other attached packages:
## [1] graphite_1.44.0 DT_0.26 cowplot_1.1.1 magrittr_2.0.3
## [5] forcats_0.5.2 stringr_1.4.1 dplyr_1.0.10 purrr_0.3.5
## [9] readr_2.1.3 tidyr_1.2.1 tibble_3.1.8 ggplot2_3.4.0
## [13] tidyverse_1.3.2 sSNAPPY_1.0.2 BiocStyle_2.26.0
##
## loaded via a namespace (and not attached):
## [1] readxl_1.4.1 backports_1.4.1
## [3] systemfonts_1.0.4 plyr_1.8.8
## [5] igraph_1.3.5 splines_4.2.2
## [7] BiocParallel_1.32.1 crosstalk_1.2.0
## [9] GenomeInfoDb_1.34.3 digest_0.6.30
## [11] htmltools_0.5.3 viridis_0.6.2
## [13] fansi_1.0.3 memoise_2.0.1
## [15] googlesheets4_1.0.1 tzdb_0.3.0
## [17] limma_3.54.0 Biostrings_2.66.0
## [19] graphlayouts_0.8.3 modelr_0.1.10
## [21] matrixStats_0.63.0 timechange_0.1.1
## [23] pkgdown_2.0.6 colorspace_2.0-3
## [25] blob_1.2.3 rvest_1.0.3
## [27] rappdirs_0.3.3 ggrepel_0.9.2
## [29] textshaping_0.3.6 haven_2.5.1
## [31] xfun_0.35 crayon_1.5.2
## [33] RCurl_1.98-1.9 jsonlite_1.8.3
## [35] graph_1.76.0 glue_1.6.2
## [37] polyclip_1.10-4 gtable_0.3.1
## [39] gargle_1.2.1 zlibbioc_1.44.0
## [41] XVector_0.38.0 DelayedArray_0.24.0
## [43] BiocGenerics_0.44.0 scales_1.2.1
## [45] pheatmap_1.0.12 DBI_1.1.3
## [47] edgeR_3.40.0 Rcpp_1.0.9
## [49] viridisLite_0.4.1 bit_4.0.5
## [51] stats4_4.2.2 htmlwidgets_1.5.4
## [53] httr_1.4.4 RColorBrewer_1.1-3
## [55] ellipsis_0.3.2 pkgconfig_2.0.3
## [57] farver_2.1.1 sass_0.4.3
## [59] dbplyr_2.2.1 locfit_1.5-9.6
## [61] utf8_1.2.2 labeling_0.4.2
## [63] tidyselect_1.2.0 rlang_1.0.6
## [65] reshape2_1.4.4 AnnotationDbi_1.60.0
## [67] munsell_0.5.0 cellranger_1.1.0
## [69] tools_4.2.2 cachem_1.0.6
## [71] cli_3.4.1 generics_0.1.3
## [73] RSQLite_2.2.18 broom_1.0.1
## [75] evaluate_0.18 fastmap_1.1.0
## [77] yaml_2.3.6 ragg_1.2.4
## [79] org.Hs.eg.db_3.16.0 knitr_1.41
## [81] bit64_4.0.5 fs_1.5.2
## [83] tidygraph_1.2.2 pander_0.6.5
## [85] KEGGREST_1.38.0 ggraph_2.1.0
## [87] nlme_3.1-160 xml2_1.3.3
## [89] compiler_4.2.2 rstudioapi_0.14
## [91] curl_4.3.3 png_0.1-7
## [93] reprex_2.0.2 tweenr_2.0.2
## [95] bslib_0.4.1 stringi_1.7.8
## [97] highr_0.9 desc_1.4.2
## [99] lattice_0.20-45 Matrix_1.5-3
## [101] vctrs_0.5.1 pillar_1.8.1
## [103] lifecycle_1.0.3 BiocManager_1.30.19
## [105] jquerylib_0.1.4 bitops_1.0-7
## [107] GenomicRanges_1.50.1 R6_2.5.1
## [109] bookdown_0.30 gridExtra_2.3
## [111] IRanges_2.32.0 codetools_0.2-18
## [113] MASS_7.3-58.1 assertthat_0.2.1
## [115] SummarizedExperiment_1.28.0 rprojroot_2.0.3
## [117] withr_2.5.0 S4Vectors_0.36.0
## [119] GenomeInfoDbData_1.2.9 mgcv_1.8-41
## [121] parallel_4.2.2 hms_1.1.2
## [123] grid_4.2.2 rmarkdown_2.18
## [125] MatrixGenerics_1.10.0 googledrive_2.0.0
## [127] ggnewscale_0.4.8 ggforce_0.4.1
## [129] Biobase_2.58.0 lubridate_1.9.0