vignettes/DeconvoBuddies.Rmd
DeconvoBuddies.Rmd
DeconvoBuddies
R
is an open-source statistical environment which can be
easily modified to enhance its functionality via packages. DeconvoBuddies
is a R
package available via the Bioconductor repository for packages.
R
can be installed on any operating system from CRAN after which you can install
DeconvoBuddies
by using the following commands in your R
session:
if (!requireNamespace("BiocManager", quietly = TRUE)) {
install.packages("BiocManager")
}
BiocManager::install("DeconvoBuddies")
## Check that you have a valid Bioconductor installation
BiocManager::valid()
DeconvoBuddies is based on many other packages and in particular in those that have implemented the infrastructure needed for dealing with snRNA-seq data. That is, packages like SingleCellExperiment.
If you are asking yourself the question “Where do I start using Bioconductor?” you might be interested in this blog post.
As package developers, we try to explain clearly how to use our
packages and in which order to use the functions. But R
and
Bioconductor
have a steep learning curve so it is critical
to learn where to ask for help. The blog post quoted above mentions some
but we would like to highlight the Bioconductor support site
as the main resource for getting help: remember to use the
DeconvoBuddies
tag and check the older
posts. Other alternatives are available such as creating GitHub
issues and tweeting. However, please note that if you want to receive
help you should adhere to the posting
guidelines. It is particularly critical that you provide a small
reproducible example and your session information so package developers
can track down the source of the error.
DeconvoBuddies
We hope that DeconvoBuddies will be useful for your research. Please use the following information to cite the package and the overall approach. Thank you!
## Citation info
citation("DeconvoBuddies")
#> To cite package 'DeconvoBuddies' in publications use:
#>
#> Huuki-Myers LA, Maynard KR, Hicks SC, Zandi P, Kleinman JE, Hyde TM,
#> Goes FS, Collado-Torres L (2024). _DeconvoBuddies: a R/Bioconductor
#> package with deconvolution helper functions_.
#> doi:10.18129/B9.bioc.DeconvoBuddies
#> <https://doi.org/10.18129/B9.bioc.DeconvoBuddies>,
#> https://github.com/LieberInstitute/DeconvoBuddies/DeconvoBuddies - R
#> package version 0.99.8,
#> <http://www.bioconductor.org/packages/DeconvoBuddies>.
#>
#> Huuki-Myers LA, Montgomery KD, Kwon SH, Cinquemani S, Eagles NJ,
#> Gonzalez-Padilla D, Maden SK, Kleinman JE, Hyde TM, Hicks SC, Maynard
#> KR, Collado-Torres L (2024). "Benchmark of cellular deconvolution
#> methods using a multi-assay reference dataset from postmortem human
#> prefrontal cortex." _bioRxiv_. doi:10.1101/2024.02.09.579665
#> <https://doi.org/10.1101/2024.02.09.579665>,
#> <https://doi.org/10.1101/2024.02.09.579665>.
#>
#> To see these entries in BibTeX format, use 'print(<citation>,
#> bibtex=TRUE)', 'toBibtex(.)', or set
#> 'options(citation.bibtex.max=999)'.
DeconvoBuddies
Let’s load some packages we’ll use in this vignette.
suppressMessages({
library("DeconvoBuddies")
library("SummarizedExperiment")
library("dplyr")
library("tidyr")
library("tibble")
})
Use fetch_deconvo_data
to download RNA sequencing data
from the Human
DLPFC (Huuki-Myers, Montgomery, Kwon, Cinquemani, Eagles,
Gonzalez-Padilla, Maden, Kleinman, Hyde, Hicks, Maynard, and
Collado-Torres, 2024).
rse_gene
: 110 samples of bulk RNA-seq. [110 bulk
RNA-seq samples x 21k genes] (41 MB).
sce
: snRNA-seq data from the Human DLPFC. [77k
nuclei x 36k genes] (172 MB)
sce_DLPFC_example
: Sub-set of sce
useful for testing. [10k nuclei x 557 genes] (49 MB)
## Access and snRNA-seq example data
if (!exists("sce_DLPFC_example")) sce_DLPFC_example <- fetch_deconvo_data("sce_DLPFC_example")
#> 2024-10-31 13:27:50.301286 loading file /github/home/.cache/R/BiocFileCache/66141272054e_sce_DLPFC_example.Rdata%3Frlkey%3Dv3z4u8ru0d2y12zgdl1az07q9%26st%3D1dcfqc1i%26dl%3D1
## Explore snRNA-seq data in sce_DLPFC_example
sce_DLPFC_example
#> class: SingleCellExperiment
#> dim: 557 10000
#> metadata(3): Samples cell_type_colors cell_type_colors_broad
#> assays(1): logcounts
#> rownames(557): GABRD PRDM16 ... AFF2 MAMLD1
#> rowData names(7): source type ... gene_type binomial_deviance
#> colnames(10000): 8_AGTGACTGTAGTTACC-1 17_GCAGCCAGTGAGTCAG-1 ...
#> 12_GGACGTCTCTGACAGT-1 1_GGTTAACTCTCTCTAA-1
#> colData names(32): Sample Barcode ... cellType_layer layer_annotation
#> reducedDimNames(0):
#> mainExpName: NULL
#> altExpNames(0):
## Access Bulk RNA-seq data
if (!exists("rse_gene")) rse_gene <- fetch_deconvo_data("rse_gene")
#> 2024-10-31 13:27:52.36966 loading file /github/home/.cache/R/BiocFileCache/66143de2ada5_rse_gene.Rdata%3Frlkey%3Dsw2djr71y954yw4o3xrmjv59b%26dl%3D1
## Explore bulk data in rse_gene
rse_gene
#> class: RangedSummarizedExperiment
#> dim: 21745 110
#> metadata(1): SPEAQeasy_settings
#> assays(2): counts logcounts
#> rownames(21745): ENSG00000227232.5 ENSG00000278267.1 ...
#> ENSG00000210195.2 ENSG00000210196.2
#> rowData names(11): Length gencodeID ... gencodeTx passExprsCut
#> colnames(110): 2107UNHS-0291_Br2720_Mid_Bulk
#> 2107UNHS-0291_Br2720_Mid_Cyto ... AN00000906_Br8667_Mid_Cyto
#> AN00000906_Br8667_Mid_Nuc
#> colData names(80): SAMPLE_ID Sample ... diagnosis qc_class
For more details on this dataset, and an example deconvolution run check out the Vignette: Deconvolution Benchmark in Human DLPFC.
Accurate deconvolution requires highly specific marker genes for each
cell type to be defined. To select genes specific for each cell type,
you can evaluate the MeanRatio for each gene x each cell type,
where
MeanRatio = mean(Expression of target cell type) / mean(Expression of highest non-target cell type)
.
These values can be calculated for a single cell RNA-seq dataset
using get_mean_ratio()
. This can also work for
spatially-resolved transcriptomics datasets. That is,
get_mean_ratio()
can also work with
SpatialExperiment::SpatialExperiment()
objects.
## find marker genes with get_mean_ratio
marker_stats <- get_mean_ratio(
sce_DLPFC_example,
cellType_col = "cellType_broad_hc",
gene_name = "gene_name",
gene_ensembl = "gene_id"
)
## explore tibble output, gene with high MeanRatio values are good marker genes
marker_stats
#> # A tibble: 762 × 10
#> gene cellType.target mean.target cellType.2nd mean.2nd MeanRatio
#> <chr> <fct> <dbl> <fct> <dbl> <dbl>
#> 1 CD22 Oligo 1.36 OPC 0.0730 18.6
#> 2 LINC01608 Oligo 2.39 Micro 0.142 16.8
#> 3 FOLH1 Oligo 1.59 OPC 0.101 15.7
#> 4 SLC5A11 Oligo 2.14 Micro 0.145 14.7
#> 5 AC012494.1 Oligo 2.42 OPC 0.169 14.3
#> 6 ST18 Oligo 4.65 OPC 0.329 14.1
#> 7 MAG Oligo 1.44 Astro 0.103 14.0
#> 8 ANLN Oligo 1.60 Micro 0.115 13.9
#> 9 CLDN11 Oligo 1.82 EndoMural 0.146 12.5
#> 10 MOG Oligo 2.06 OPC 0.185 11.1
#> # ℹ 752 more rows
#> # ℹ 4 more variables: MeanRatio.rank <int>, MeanRatio.anno <chr>,
#> # gene_ensembl <chr>, gene_name <chr>
For more discussion of finding marker genes with
DeconvoBuddies
check out the Vignette:
Finding Marker Genes with DeconvoBuddies.
As you work with single-cell data and deconvolution outputs, it is
very useful to establish a consistent color pallet to use across
different plots. The function create_cell_colors()
returns
a named vector of hex values, corresponding to the names of cell types.
This list is compatible with functions like
ggplot2::scale_color_manual()
.
There are three pallets to choose from to generate colors:
“classic” (default): "Set1"
from
RColorBrewer
- max 9 colors
“gg”: Equi-distant hues, same process for selecting colors as
ggplot
- no maximum number
“tableau”: tableau20 color set - max 20 colors
test_cell_types <- c("cell_A", "cell_B", "cell_C", "cell_D", "cell_E")
## Preview "classic" colors
test_cell_colors_classic <- create_cell_colors(
cell_types = test_cell_types,
pallet = "classic",
preview = TRUE
)
## Preview "gg" colors
test_cell_colors_gg <- create_cell_colors(
cell_types = test_cell_types,
pallet = "gg",
preview = TRUE
)
## Preview "tableau" colors
test_cell_colors_tableau <- create_cell_colors(
cell_types = test_cell_types,
pallet = "tableau",
preview = TRUE
)
## Check the color hex codes for "tableau"
test_cell_colors_tableau
#> cell_A cell_B cell_C cell_D cell_E
#> "#1F77B4" "#AEC7E8" "#FF7F0E" "#FFBB78" "#2CA02C"
If there are sub-cell types with consistent delimiters, the
split
argument creates a scale of related colors. This
helps expand on the maximum number of colors and makes your pallet
flexible when considering different ‘resolutions’ of cell types. This
works by ignoring any prefixes after the split
character.
In this example below, Excit_01
and Excit_02
will just be considered as Excit
since
split = "_"
.
my_cell_types <- levels(sce_DLPFC_example$cellType_hc)
## Ignore any suffix after the "_" character by using the "split" argument
my_cell_colors <- create_cell_colors(
cell_types = my_cell_types,
pallet = "classic",
preview = TRUE,
split = "_"
)
The function plot_marker_express()
helps quickly
visualize expression of top marker genes, by ordering and annotating
violin plots of expression over cell type. Here we’ll plot the
expression of the top 6 marker genes for Astrocytes.
# plot expression of the top 6 Astro marker genes
plot_marker_express(
sce = sce_DLPFC_example,
stats = marker_stats,
cell_type = "Astro",
n_genes = 6,
cellType_col = "cellType_broad_hc",
color_pal = my_cell_colors
)
The violin plots of gene expression confirm the cell type specificity of these marker genes, most of the nuclei with high expression of these six genes are astrocytes (Astro).
The output of deconvolution are cell type estimates that sum to 1. A
good visulization for these predictions is a stacked bar plot. The
function plot_composition_bar()
creates a stacked bar plot
showing the cell type proportion for each sample, or the average
proportion for a group of samples. In this example data, the
RNum
is a sample (donor) identifier and Dx
is
a group variable for the diagnosis status of the donors.
# load example data
data("rse_bulk_test")
data("est_prop")
# access the colData of a test rse dataset
pd <- colData(rse_bulk_test) |>
as.data.frame()
## pivot data to long format and join with test estimated proportion data
est_prop_long <- est_prop |>
rownames_to_column("RNum") |>
pivot_longer(!RNum, names_to = "cell_type", values_to = "prop") |>
left_join(pd)
#> Joining with `by = join_by(RNum)`
## explore est_prop_long
est_prop_long
#> # A tibble: 500 × 7
#> RNum cell_type prop BrNum Sex Dx Age
#> <chr> <chr> <dbl> <chr> <chr> <chr> <dbl>
#> 1 R913 cell_A 0.380 Br001 F Case 71.9
#> 2 R913 cell_B 0.240 Br001 F Case 71.9
#> 3 R913 cell_C 0.162 Br001 F Case 71.9
#> 4 R913 cell_D 0.108 Br001 F Case 71.9
#> 5 R913 cell_E 0.111 Br001 F Case 71.9
#> 6 R602 cell_A 0.443 Br002 F Control 73.1
#> 7 R602 cell_B 0.207 Br002 F Control 73.1
#> 8 R602 cell_C 0.0743 Br002 F Control 73.1
#> 9 R602 cell_D 0.160 Br002 F Control 73.1
#> 10 R602 cell_E 0.115 Br002 F Control 73.1
#> # ℹ 490 more rows
## the composition bar plot shows cell type composition for Sample
plot_composition_bar(est_prop_long,
x_col = "RNum",
add_text = FALSE
) +
ggplot2::scale_fill_manual(values = test_cell_colors_classic)
## the composition bar plot shows the average cell type composition for each Dx
plot_composition_bar(est_prop_long, x_col = "Dx") +
ggplot2::scale_fill_manual(values = test_cell_colors_classic)
We can see that the mean proportions of cell types A through E are
very similar across the Dx
groups (Case
and
Control
). In this case, this is expected given that we are
using simulated data. Although if you look across each donor with
RNum
we can see more variability across the simulated
data.
Since you are now familiar with the basic overview of
DeconvoBuddies
, you are now ready to dive deeper into:
DeconvoBuddies
to the Human
Brain (DLPFC) Deconvolution dataset.The DeconvoBuddies package (Huuki-Myers, Maynard, Hicks, Zandi, Kleinman, Hyde, Goes, and Collado-Torres, 2024) was made possible thanks to:
This package was developed using biocthis.
Code for creating the vignette
## Create the vignette
library("rmarkdown")
system.time(render("DeconvoBuddies.Rmd", "BiocStyle::html_document"))
## Extract the R code
library("knitr")
knit("DeconvoBuddies.Rmd", tangle = TRUE)
Date the vignette was generated.
#> [1] "2024-10-31 13:27:57 UTC"
Wallclock time spent generating the vignette.
#> Time difference of 19.833 secs
R
session information.
#> ─ Session info ───────────────────────────────────────────────────────────────────────────────────────────────────────
#> setting value
#> version R version 4.4.1 (2024-06-14)
#> os Ubuntu 22.04.5 LTS
#> system x86_64, linux-gnu
#> ui X11
#> language en
#> collate en_US.UTF-8
#> ctype en_US.UTF-8
#> tz UTC
#> date 2024-10-31
#> pandoc 3.4 @ /usr/bin/ (via rmarkdown)
#>
#> ─ Packages ───────────────────────────────────────────────────────────────────────────────────────────────────────────
#> package * version date (UTC) lib source
#> abind 1.4-8 2024-09-12 [1] RSPM (R 4.4.0)
#> AnnotationDbi 1.68.0 2024-10-29 [1] Bioconductor 3.20 (R 4.4.1)
#> AnnotationHub 3.14.0 2024-10-29 [1] Bioconductor 3.20 (R 4.4.1)
#> attempt 0.3.1 2020-05-03 [1] RSPM (R 4.4.0)
#> backports 1.5.0 2024-05-23 [1] RSPM (R 4.4.0)
#> beachmat 2.22.0 2024-10-29 [1] Bioconductor 3.20 (R 4.4.1)
#> beeswarm 0.4.0 2021-06-01 [1] RSPM (R 4.4.0)
#> benchmarkme 1.0.8 2022-06-12 [1] RSPM (R 4.4.0)
#> benchmarkmeData 1.0.4 2020-04-23 [1] RSPM (R 4.4.0)
#> bibtex 0.5.1 2023-01-26 [1] RSPM (R 4.4.0)
#> Biobase * 2.66.0 2024-10-29 [1] Bioconductor 3.20 (R 4.4.1)
#> BiocFileCache 2.14.0 2024-10-29 [1] Bioconductor 3.20 (R 4.4.1)
#> BiocGenerics * 0.52.0 2024-10-29 [1] Bioconductor 3.20 (R 4.4.1)
#> BiocIO 1.16.0 2024-10-29 [1] Bioconductor 3.20 (R 4.4.1)
#> BiocManager 1.30.25 2024-08-28 [2] CRAN (R 4.4.1)
#> BiocNeighbors 2.0.0 2024-10-29 [1] Bioconductor 3.20 (R 4.4.1)
#> BiocParallel 1.40.0 2024-10-29 [1] Bioconductor 3.20 (R 4.4.1)
#> BiocSingular 1.22.0 2024-10-29 [1] Bioconductor 3.20 (R 4.4.1)
#> BiocStyle * 2.34.0 2024-10-29 [1] Bioconductor 3.20 (R 4.4.1)
#> BiocVersion 3.20.0 2024-10-21 [2] Bioconductor 3.20 (R 4.4.1)
#> Biostrings 2.74.0 2024-10-29 [1] Bioconductor 3.20 (R 4.4.1)
#> bit 4.5.0 2024-09-20 [1] RSPM (R 4.4.0)
#> bit64 4.5.2 2024-09-22 [1] RSPM (R 4.4.0)
#> bitops 1.0-9 2024-10-03 [1] RSPM (R 4.4.0)
#> blob 1.2.4 2023-03-17 [1] RSPM (R 4.4.0)
#> bluster 1.16.0 2024-10-29 [1] Bioconductor 3.20 (R 4.4.1)
#> bookdown 0.41 2024-10-16 [1] RSPM (R 4.4.0)
#> bslib 0.8.0 2024-07-29 [2] RSPM (R 4.4.0)
#> cachem 1.1.0 2024-05-16 [2] RSPM (R 4.4.0)
#> cli 3.6.3 2024-06-21 [2] RSPM (R 4.4.0)
#> cluster 2.1.6 2023-12-01 [3] CRAN (R 4.4.1)
#> codetools 0.2-20 2024-03-31 [3] CRAN (R 4.4.1)
#> colorspace 2.1-1 2024-07-26 [1] RSPM (R 4.4.0)
#> config 0.3.2 2023-08-30 [1] RSPM (R 4.4.0)
#> cowplot 1.1.3 2024-01-22 [1] RSPM (R 4.4.0)
#> crayon 1.5.3 2024-06-20 [2] RSPM (R 4.4.0)
#> curl 5.2.3 2024-09-20 [2] RSPM (R 4.4.0)
#> data.table 1.16.2 2024-10-10 [1] RSPM (R 4.4.0)
#> DBI 1.2.3 2024-06-02 [1] RSPM (R 4.4.0)
#> dbplyr 2.5.0 2024-03-19 [1] RSPM (R 4.4.0)
#> DeconvoBuddies * 0.99.8 2024-10-31 [1] Bioconductor
#> DelayedArray 0.32.0 2024-10-29 [1] Bioconductor 3.20 (R 4.4.1)
#> desc 1.4.3 2023-12-10 [2] RSPM (R 4.4.0)
#> digest 0.6.37 2024-08-19 [2] RSPM (R 4.4.0)
#> doParallel 1.0.17 2022-02-07 [1] RSPM (R 4.4.0)
#> dotCall64 1.2 2024-10-04 [1] RSPM (R 4.4.0)
#> dplyr * 1.1.4 2023-11-17 [1] RSPM (R 4.4.0)
#> dqrng 0.4.1 2024-05-28 [1] RSPM (R 4.4.0)
#> DT 0.33 2024-04-04 [1] RSPM (R 4.4.0)
#> edgeR 4.4.0 2024-10-29 [1] Bioconductor 3.20 (R 4.4.1)
#> evaluate 1.0.1 2024-10-10 [2] RSPM (R 4.4.0)
#> ExperimentHub 2.14.0 2024-10-29 [1] Bioconductor 3.20 (R 4.4.1)
#> fansi 1.0.6 2023-12-08 [2] RSPM (R 4.4.0)
#> farver 2.1.2 2024-05-13 [1] RSPM (R 4.4.0)
#> fastmap 1.2.0 2024-05-15 [2] RSPM (R 4.4.0)
#> fields 16.3 2024-09-30 [1] RSPM (R 4.4.0)
#> filelock 1.0.3 2023-12-11 [1] RSPM (R 4.4.0)
#> foreach 1.5.2 2022-02-02 [1] RSPM (R 4.4.0)
#> fs 1.6.5 2024-10-30 [2] RSPM (R 4.4.0)
#> generics 0.1.3 2022-07-05 [1] RSPM (R 4.4.0)
#> GenomeInfoDb * 1.42.0 2024-10-29 [1] Bioconductor 3.20 (R 4.4.1)
#> GenomeInfoDbData 1.2.13 2024-10-11 [1] Bioconductor
#> GenomicAlignments 1.42.0 2024-10-29 [1] Bioconductor 3.20 (R 4.4.1)
#> GenomicRanges * 1.58.0 2024-10-29 [1] Bioconductor 3.20 (R 4.4.1)
#> ggbeeswarm 0.7.2 2023-04-29 [1] RSPM (R 4.4.0)
#> ggplot2 3.5.1 2024-04-23 [1] RSPM (R 4.4.0)
#> ggrepel 0.9.6 2024-09-07 [1] RSPM (R 4.4.0)
#> glue 1.8.0 2024-09-30 [2] RSPM (R 4.4.0)
#> golem 0.5.1 2024-08-27 [1] RSPM (R 4.4.0)
#> gridExtra 2.3 2017-09-09 [1] RSPM (R 4.4.0)
#> gtable 0.3.6 2024-10-25 [1] RSPM (R 4.4.0)
#> highr 0.11 2024-05-26 [2] RSPM (R 4.4.0)
#> htmltools 0.5.8.1 2024-04-04 [2] RSPM (R 4.4.0)
#> htmlwidgets 1.6.4 2023-12-06 [2] RSPM (R 4.4.0)
#> httpuv 1.6.15 2024-03-26 [2] RSPM (R 4.4.0)
#> httr 1.4.7 2023-08-15 [1] RSPM (R 4.4.0)
#> igraph 2.1.1 2024-10-19 [1] RSPM (R 4.4.0)
#> IRanges * 2.40.0 2024-10-29 [1] Bioconductor 3.20 (R 4.4.1)
#> irlba 2.3.5.1 2022-10-03 [1] RSPM (R 4.4.0)
#> iterators 1.0.14 2022-02-05 [1] RSPM (R 4.4.0)
#> jquerylib 0.1.4 2021-04-26 [2] RSPM (R 4.4.0)
#> jsonlite 1.8.9 2024-09-20 [2] RSPM (R 4.4.0)
#> KEGGREST 1.46.0 2024-10-29 [1] Bioconductor 3.20 (R 4.4.1)
#> knitr 1.48 2024-07-07 [2] RSPM (R 4.4.0)
#> labeling 0.4.3 2023-08-29 [1] RSPM (R 4.4.0)
#> later 1.3.2 2023-12-06 [2] RSPM (R 4.4.0)
#> lattice 0.22-6 2024-03-20 [3] CRAN (R 4.4.1)
#> lazyeval 0.2.2 2019-03-15 [1] RSPM (R 4.4.0)
#> lifecycle 1.0.4 2023-11-07 [2] RSPM (R 4.4.0)
#> limma 3.62.0 2024-10-29 [1] Bioconductor 3.20 (R 4.4.1)
#> locfit 1.5-9.10 2024-06-24 [1] RSPM (R 4.4.0)
#> lubridate 1.9.3 2023-09-27 [1] RSPM (R 4.4.0)
#> magick 2.8.5 2024-09-20 [1] RSPM (R 4.4.0)
#> magrittr 2.0.3 2022-03-30 [2] RSPM (R 4.4.0)
#> maps 3.4.2 2023-12-15 [1] RSPM (R 4.4.0)
#> Matrix 1.7-1 2024-10-18 [2] RSPM (R 4.4.0)
#> MatrixGenerics * 1.18.0 2024-10-29 [1] Bioconductor 3.20 (R 4.4.1)
#> matrixStats * 1.4.1 2024-09-08 [1] RSPM (R 4.4.0)
#> memoise 2.0.1 2021-11-26 [2] RSPM (R 4.4.0)
#> metapod 1.14.0 2024-10-29 [1] Bioconductor 3.20 (R 4.4.1)
#> mime 0.12 2021-09-28 [2] RSPM (R 4.4.0)
#> munsell 0.5.1 2024-04-01 [1] RSPM (R 4.4.0)
#> paletteer 1.6.0 2024-01-21 [1] RSPM (R 4.4.0)
#> pillar 1.9.0 2023-03-22 [2] RSPM (R 4.4.0)
#> pkgconfig 2.0.3 2019-09-22 [2] RSPM (R 4.4.0)
#> pkgdown 2.1.1 2024-09-17 [1] RSPM (R 4.4.0)
#> plotly 4.10.4 2024-01-13 [1] RSPM (R 4.4.0)
#> plyr 1.8.9 2023-10-02 [1] RSPM (R 4.4.0)
#> png 0.1-8 2022-11-29 [1] RSPM (R 4.4.0)
#> promises 1.3.0 2024-04-05 [2] RSPM (R 4.4.0)
#> purrr 1.0.2 2023-08-10 [2] RSPM (R 4.4.0)
#> R6 2.5.1 2021-08-19 [2] RSPM (R 4.4.0)
#> rafalib 1.0.0 2015-08-09 [1] RSPM (R 4.4.0)
#> ragg 1.3.3 2024-09-11 [2] RSPM (R 4.4.0)
#> rappdirs 0.3.3 2021-01-31 [2] RSPM (R 4.4.0)
#> RColorBrewer 1.1-3 2022-04-03 [1] RSPM (R 4.4.0)
#> Rcpp 1.0.13 2024-07-17 [2] RSPM (R 4.4.0)
#> RCurl 1.98-1.16 2024-07-11 [1] RSPM (R 4.4.0)
#> RefManageR * 1.4.0 2022-09-30 [1] RSPM (R 4.4.0)
#> rematch2 2.1.2 2020-05-01 [2] RSPM (R 4.4.0)
#> reshape2 1.4.4 2020-04-09 [1] RSPM (R 4.4.0)
#> restfulr 0.0.15 2022-06-16 [1] RSPM (R 4.4.0)
#> rjson 0.2.23 2024-09-16 [1] RSPM (R 4.4.0)
#> rlang 1.1.4 2024-06-04 [2] RSPM (R 4.4.0)
#> rmarkdown 2.28 2024-08-17 [2] RSPM (R 4.4.0)
#> Rsamtools 2.22.0 2024-10-29 [1] Bioconductor 3.20 (R 4.4.1)
#> RSQLite 2.3.7 2024-05-27 [1] RSPM (R 4.4.0)
#> rsvd 1.0.5 2021-04-16 [1] RSPM (R 4.4.0)
#> rtracklayer 1.66.0 2024-10-29 [1] Bioconductor 3.20 (R 4.4.1)
#> S4Arrays 1.6.0 2024-10-29 [1] Bioconductor 3.20 (R 4.4.1)
#> S4Vectors * 0.44.0 2024-10-29 [1] Bioconductor 3.20 (R 4.4.1)
#> sass 0.4.9 2024-03-15 [2] RSPM (R 4.4.0)
#> ScaledMatrix 1.14.0 2024-10-29 [1] Bioconductor 3.20 (R 4.4.1)
#> scales 1.3.0 2023-11-28 [1] RSPM (R 4.4.0)
#> scater 1.34.0 2024-10-29 [1] Bioconductor 3.20 (R 4.4.1)
#> scran 1.34.0 2024-10-29 [1] Bioconductor 3.20 (R 4.4.1)
#> scuttle 1.16.0 2024-10-29 [1] Bioconductor 3.20 (R 4.4.1)
#> sessioninfo * 1.2.2 2021-12-06 [2] RSPM (R 4.4.0)
#> shiny 1.9.1 2024-08-01 [2] RSPM (R 4.4.0)
#> shinyWidgets 0.8.7 2024-09-23 [1] RSPM (R 4.4.0)
#> SingleCellExperiment 1.28.0 2024-10-29 [1] Bioconductor 3.20 (R 4.4.1)
#> spam 2.11-0 2024-10-03 [1] RSPM (R 4.4.0)
#> SparseArray 1.6.0 2024-10-29 [1] Bioconductor 3.20 (R 4.4.1)
#> SpatialExperiment 1.16.0 2024-10-29 [1] Bioconductor 3.20 (R 4.4.1)
#> spatialLIBD 1.17.8 2024-07-25 [1] Bioconductor 3.20 (R 4.4.1)
#> statmod 1.5.0 2023-01-06 [1] RSPM (R 4.4.0)
#> stringi 1.8.4 2024-05-06 [2] RSPM (R 4.4.0)
#> stringr 1.5.1 2023-11-14 [2] RSPM (R 4.4.0)
#> SummarizedExperiment * 1.36.0 2024-10-29 [1] Bioconductor 3.20 (R 4.4.1)
#> systemfonts 1.1.0 2024-05-15 [2] RSPM (R 4.4.0)
#> textshaping 0.4.0 2024-05-24 [2] RSPM (R 4.4.0)
#> tibble * 3.2.1 2023-03-20 [2] RSPM (R 4.4.0)
#> tidyr * 1.3.1 2024-01-24 [1] RSPM (R 4.4.0)
#> tidyselect 1.2.1 2024-03-11 [1] RSPM (R 4.4.0)
#> timechange 0.3.0 2024-01-18 [1] RSPM (R 4.4.0)
#> UCSC.utils 1.2.0 2024-10-29 [1] Bioconductor 3.20 (R 4.4.1)
#> utf8 1.2.4 2023-10-22 [2] RSPM (R 4.4.0)
#> vctrs 0.6.5 2023-12-01 [2] RSPM (R 4.4.0)
#> vipor 0.4.7 2023-12-18 [1] RSPM (R 4.4.0)
#> viridis 0.6.5 2024-01-29 [1] RSPM (R 4.4.0)
#> viridisLite 0.4.2 2023-05-02 [1] RSPM (R 4.4.0)
#> withr 3.0.2 2024-10-28 [2] RSPM (R 4.4.0)
#> xfun 0.48 2024-10-03 [2] RSPM (R 4.4.0)
#> XML 3.99-0.17 2024-06-25 [1] RSPM (R 4.4.0)
#> xml2 1.3.6 2023-12-04 [2] RSPM (R 4.4.0)
#> xtable 1.8-4 2019-04-21 [2] RSPM (R 4.4.0)
#> XVector 0.46.0 2024-10-29 [1] Bioconductor 3.20 (R 4.4.1)
#> yaml 2.3.10 2024-07-26 [2] RSPM (R 4.4.0)
#> zlibbioc 1.52.0 2024-10-29 [1] Bioconductor 3.20 (R 4.4.1)
#>
#> [1] /__w/_temp/Library
#> [2] /usr/local/lib/R/site-library
#> [3] /usr/local/lib/R/library
#>
#> ──────────────────────────────────────────────────────────────────────────────────────────────────────────────────────
This vignette was generated using BiocStyle (Oleś, 2024) with knitr (Xie, 2024) and rmarkdown (Allaire, Xie, Dervieux et al., 2024) running behind the scenes.
Citations made with RefManageR (McLean, 2017).
[1] J. Allaire, Y. Xie, C. Dervieux, et al. rmarkdown: Dynamic Documents for R. R package version 2.28. 2024. URL: https://github.com/rstudio/rmarkdown.
[2] L. A. Huuki-Myers, K. R. Maynard, S. C. Hicks, et al. DeconvoBuddies: a R/Bioconductor package with deconvolution helper functions. https://github.com/LieberInstitute/DeconvoBuddies/DeconvoBuddies - R package version 0.99.8. 2024. DOI: 10.18129/B9.bioc.DeconvoBuddies. URL: http://www.bioconductor.org/packages/DeconvoBuddies.
[3] L. A. Huuki-Myers, K. D. Montgomery, S. H. Kwon, et al. “Benchmark of cellular deconvolution methods using a multi-assay reference dataset from postmortem human prefrontal cortex”. In: bioRxiv (2024). DOI: 10.1101/2024.02.09.579665. URL: https://doi.org/10.1101/2024.02.09.579665.
[4] M. W. McLean. “RefManageR: Import and Manage BibTeX and BibLaTeX References in R”. In: The Journal of Open Source Software (2017). DOI: 10.21105/joss.00338.
[5] A. Oleś. BiocStyle: Standard styles for vignettes and other Bioconductor documents. R package version 2.34.0. 2024. DOI: 10.18129/B9.bioc.BiocStyle. URL: https://bioconductor.org/packages/BiocStyle.
[6] R Core Team. R: A Language and Environment for Statistical Computing. R Foundation for Statistical Computing. Vienna, Austria, 2024. URL: https://www.R-project.org/.
[7] H. Wickham. “testthat: Get Started with Testing”. In: The R Journal 3 (2011), pp. 5–10. URL: https://journal.r-project.org/archive/2011-1/RJournal_2011-1_Wickham.pdf.
[8] H. Wickham, W. Chang, R. Flight, et al. sessioninfo: R Session Information. R package version 1.2.2, https://r-lib.github.io/sessioninfo/. 2021. URL: https://github.com/r-lib/sessioninfo#readme.
[9] Y. Xie. knitr: A General-Purpose Package for Dynamic Report Generation in R. R package version 1.48. 2024. URL: https://yihui.org/knitr/.