TREG 1.2.0
Note: TREG is pronounced as a single word and fully capitalized, unlike Regulatory T cells, which are known as “Tregs” (pronounced “T-regs”). The work described here is unrelated to regulatory T cells.
TREGR is an open-source statistical environment which can be easily modified to enhance its functionality via packages. TREG is a R package available via Bioconductor. R can be installed on any operating system from CRAN after which you can install TREG by using the following commands in your R session:
if (!requireNamespace("BiocManager", quietly = TRUE)) {
    install.packages("BiocManager")
}
BiocManager::install("TREG")
## Check that you have a valid Bioconductor installation
BiocManager::valid()TREG (Huuki-Myers and Collado-Torres, 2022) is based on many other packages and in particular in those that have implemented the infrastructure needed for dealing with single cell RNA sequencing data, visualization functions, and interactive data exploration. That is, packages like SummarizedExperiment that allow you to store the data.
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 regarding Bioconductor. 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.
TREGWe hope that TREG will be useful for your research. Please use the following information to cite the package and the research article describing the data provided by TREG. Thank you!
## Citation info
citation("TREG")
#> 
#> To cite package 'TREG' in publications use:
#> 
#>   Huuki-Myers LA, Collado-Torres L (2022). _TREG: a R/Bioconductor
#>   package to identify Total RNA Expression Genes_.
#>   doi:10.18129/B9.bioc.TREG <https://doi.org/10.18129/B9.bioc.TREG>,
#>   https://github.com/LieberInstitute/TREG/TREG - R package version
#>   1.2.0, <http://www.bioconductor.org/packages/TREG>.
#> 
#>   Huuki-Myers LA, Montgomery KD, Kwon SH, Page SC, Hicks SC, Maynard
#>   KR, Collado-Torres L (2022). "Data Driven Identification of Total RNA
#>   Expression Genes "TREGs" for estimation of RNA abundance in
#>   heterogeneous cell types." _bioRxiv_. doi:10.1101/2022.04.28.489923
#>   <https://doi.org/10.1101/2022.04.28.489923>,
#>   <https://doi.org/10.1101/2022.04.28.489923>.
#> 
#> To see these entries in BibTeX format, use 'print(<citation>,
#> bibtex=TRUE)', 'toBibtex(.)', or set
#> 'options(citation.bibtex.max=999)'.The TREG (Huuki-Myers and Collado-Torres, 2022) package was developed for identifying candidate Total RNA Expression Genes (TREGs) for estimating RNA abundance for individual cells in an snFISH experiment by researchers at the Lieber Institute for Brain Development (LIBD) (Huuki-Myers, Montgomery, Kwon, Page, Hicks, Maynard, and Collado-Torres, 2022).
In this vignette we’ll showcase how you can use the R functions provided by TREG (Huuki-Myers and Collado-Torres, 2022) with the snRNA-seq dataset that was recently published by our LIBD collaborators (Tran, Maynard, Spangler, Huuki, Montgomery, Sadashivaiah, Tippani, Barry, Hancock, Hicks, Kleinman, Hyde, Collado-Torres, Jaffe, and Martinowich, 2021).
To get started, please load the TREG package.
library("TREG")The goal of TREG is to help find candidate Total RNA Expression Genes (TREGs)
in single nucleus (or single cell) RNA-seq data.
The expression of a TREG is proportional to the the overall RNA expression in a cell. This relationship can be used to estimate total RNA content in cells in assays where only a few genes can be measured, such as single-molecule fluorescent in situ hybridization (smFISH).
In a smFISH experiment the number of TREG puncta can be used to infer the total RNA expression of the cell.
The motivation of this work is to collect data via smFISH in to help build better deconvolution algorithms. But may be many other application for TREGs in experimental design!
The gene must have non-zero expression in most cells across different tissue and cell types.
A TREG should also be expressed at a constant level in respect to other genes across different cell types or have high rank invariance.
Be measurable as a continuous metric in the experimental assay, for example have a dynamic range of puncta when observed in RNAscope. This will need to be considered for the candidate TREGs, and may need to be validated experimentally.
TREG
Filter for low Proportion Zero genes snRNA-seq dataset: This is
facilitated with the functions get_prop_zero() and filter_prop_zero() (Equation (1)).
snRNA-seq data is notoriously sparse, these functions enrich for genes with more
universal expression.
Evaluate genes for Rank Invariance The nuclei are grouped only
by cell type. Within each cell type, the mean expression for each
gene is ranked, the result is a vector (length is the number of
genes), using the function rank_group(). Then the expression of each gene is
ranked for each nucleus,the result is a matrix (the number of nuclei x number
of genes), using the function rank_cells().Then the absolute difference
between the rank of each nucleus and the mean expression is found, from here
the mean of the differences for each gene is calculated, then ranked.
These steps are repeated for each group, the result is a matrix of ranks, (number of cell
types x number of genes). From here the sum of the ranks for each
gene are reversed ranked, so there is one final value for each gene,
the “Rank Invariance” The genes with the highest rank-invariance are
considered good candidates as TREGs. This is calculated with rank_invariance_express().
This full process is implemented by: rank_invariance_express().
In this example we will apply our data driven process for TREG discovery to a
snRNA-seq dataset. This process has three main steps:
1. Data prep
2. Gene filtering: dropping genes with low expression and high Proportion Zero (Equation (1))
3. Rank Invariance Calculation
library("SingleCellExperiment")
library("pheatmap")
library("dplyr")
library("ggplot2")
library("tidyr")
library("tibble")Here we download a public single nucleus RNA-seq (snRNA-seq) data from (Tran, Maynard, Spangler et al., 2021) that we’ll use as our example. This data can be accessed on github. This data is from postmortem human brain in the dorsolateral prefrontal cortex (DLPFC) region, and contains gene expression data for 11k nuclei.
We will use BiocFileCache() to cache this data. It is stored as a SingleCellExperiment
object named sce.dlpfc.tran, and takes 1.01 GB of RAM memory to load.
# Download and save a local cache of the data available at:
# https://github.com/LieberInstitute/10xPilot_snRNAseq-human#processed-data
bfc <- BiocFileCache::BiocFileCache()
url <- paste0(
    "https://libd-snrnaseq-pilot.s3.us-east-2.amazonaws.com/",
    "SCE_DLPFC-n3_tran-etal.rda"
)
local_data <- BiocFileCache::bfcrpath(url, x = bfc)
load(local_data, verbose = TRUE)
#> Loading objects:
#>   sce.dlpfc.tran| Cell Type | Acronym | 
|---|---|
| Astrocyte | Astro | 
| Excitatory Neurons | Excit | 
| Microglia | Micro | 
| Oligodendrocytes | Oligo | 
| Oligodendrocyte Progenitor Cells | OPC | 
| Inhibitory Neurons | Inhib | 
Human brain tissue consists of many types of cells, for the porpose of this demo, we will focus on the six major cell types listed in Table 1.
First we will combine all of the Excit, Inhib subtypes, as it is a finer
resolution than we want to examine, and combine rare subtypes in to one group.
If there are too few cells in a group there may not be enough data to
get good results. This new cell type classification is stored in the colData as
cellType.broad.
## Explore the dimensions and cell type annotations
dim(sce.dlpfc.tran)
#> [1] 33538 11202
table(sce.dlpfc.tran$cellType)
#> 
#>      Astro    Excit_A    Excit_B    Excit_C    Excit_D    Excit_E    Excit_F 
#>        782        529        773        524        132        187        243 
#>    Inhib_A    Inhib_B    Inhib_C    Inhib_D    Inhib_E    Inhib_F Macrophage 
#>        333        454        365        413          7          8         10 
#>      Micro      Mural      Oligo        OPC      Tcell 
#>        388         18       5455        572          9
## Use a lower resolution of cell type annotation
sce.dlpfc.tran$cellType.broad <- gsub("_[A-Z]$", "", sce.dlpfc.tran$cellType)
(cell_type_tab <- table(sce.dlpfc.tran$cellType.broad))
#> 
#>      Astro      Excit      Inhib Macrophage      Micro      Mural        OPC 
#>        782       2388       1580         10        388         18        572 
#>      Oligo      Tcell 
#>       5455          9Next, we will drop any groups with < 50 cells after merging subtypes. This excludes any very rare cell types. Now we are working with the six broad cell types we are interested in.
## Find cell types with < 50 cells
(ct_drop <- names(cell_type_tab)[cell_type_tab < 50])
#> [1] "Macrophage" "Mural"      "Tcell"
## Filter columns of sce object
sce.dlpfc.tran <- sce.dlpfc.tran[, !sce.dlpfc.tran$cellType.broad %in% ct_drop]
## Check new cell type bread down and dimension
table(sce.dlpfc.tran$cellType.broad)
#> 
#> Astro Excit Inhib Micro   OPC Oligo 
#>   782  2388  1580   388   572  5455
dim(sce.dlpfc.tran)
#> [1] 33538 11165Single Nucleus data is often very sparse (lots of zeros in the count data), this dataset is 88% sparse. We can illustrate this in the heat map of the first 1k genes and 500 cells. The heatmap is mostly blue, indicating low values (Figure 1).
## this data is 88% sparse
sum(assays(sce.dlpfc.tran)$counts == 0) / (nrow(sce.dlpfc.tran) * ncol(sce.dlpfc.tran))
#> [1] 0.8839022
## lets make a heatmap of the first 1k genes and 500 cells
count_test <- as.matrix(assays(sce.dlpfc.tran)$logcounts[seq_len(1000), seq_len(500)])
pheatmap(count_test,
    cluster_rows = FALSE,
    cluster_cols = FALSE,
    show_rownames = FALSE,
    show_colnames = FALSE
)
Figure 1: Heatmap of the snRNA-seq counts
Illustrates sparseness of unfiltered data.
Determine the median expression of genes over all rows, drop all the genes that are below this limit.
row_means <- rowMeans(assays(sce.dlpfc.tran)$logcounts)
(median_row_means <- median(row_means))
#> [1] 0.01405562
sce.dlpfc.tran <- sce.dlpfc.tran[row_means > median_row_means, ]
dim(sce.dlpfc.tran)
#> [1] 16769 11165After this filter lets check sparsity and make a heatmap of the first 1k genes and 500 cells. We are seeing more non-blue (Figure 2)!
## this data down to 77% sparse
sum(assays(sce.dlpfc.tran)$counts == 0) / (nrow(sce.dlpfc.tran) * ncol(sce.dlpfc.tran))
#> [1] 0.7713423
## replot heatmap
count_test <- as.matrix(assays(sce.dlpfc.tran)$logcounts[seq_len(1000), seq_len(500)])
pheatmap::pheatmap(count_test,
    cluster_rows = FALSE,
    cluster_cols = FALSE,
    show_rownames = FALSE,
    show_colnames = FALSE
)
Figure 2: Heatmap of the snRNA-seq counts
With top 50% filtering the data becomes less sparse.
For each group (let’s use cellType.broad) get the Proportion Zero for each gene, where Proportion Zero is defined in Equation (1) where \(c_{i,j,k,z}\) is the number of snRNA-seq counts for cell/nucleus \(z\) for gene \(i\), cell type \(j\), brain region \(k\), and \(n_{j, k}\) is the number of cells/nuclei for cell type \(j\) and brain region \(k\).
\[\begin{equation} PZ_{i,j,k} = \sum_{z=1}^{n_{j,k}} I(c_{i,j,k} > 0) / n_{j,k} \tag{1} \end{equation}\]
# get prop zero for each gene for each cell type
prop_zeros <- get_prop_zero(sce.dlpfc.tran, group_col = "cellType.broad")
head(prop_zeros)
#>                Astro     Excit     Inhib     Micro       OPC     Oligo
#> AL627309.1 0.9667519 0.7784757 0.8487342 0.9536082 0.9685315 0.9552704
#> AL669831.5 0.8836317 0.2361809 0.4879747 0.8608247 0.7657343 0.8445463
#> LINC00115  0.9948849 0.9493300 0.9759494 0.9845361 0.9930070 0.9891842
#> NOC2L      0.9347826 0.6168342 0.7556962 0.9536082 0.8811189 0.9365720
#> KLHL17     0.9859335 0.8986600 0.9360759 0.9922680 0.9877622 0.9963336
#> HES4       0.8695652 0.7357621 0.6892405 0.9974227 0.9737762 0.9906508To determine a good cutoff for filtering lets examine the distribution of these Proportion Zeros by group.
# Pivot data longer for plotting
prop_zero_long <- prop_zeros %>%
    rownames_to_column("Gene") %>%
    pivot_longer(!Gene, names_to = "Group", values_to = "prop_zero")
# Plot histograms
(prop_zero_histogram <- ggplot(
    data = prop_zero_long,
    aes(x = prop_zero, fill = Group)
) +
    geom_histogram(binwidth = 0.05) +
    facet_wrap(~Group))
Figure 3: Distribution of Proption Zero across cell types and regions
Looks like around 0.9 the densities peak, we’ll set that as the cutoff (Figure 3).
## Specify a cutoff, here we use 0.9
propZero_limit <- 0.9
## Add a vertical red dashed line where the cutoff is located
prop_zero_histogram +
    geom_vline(xintercept = propZero_limit, color = "red", linetype = "dashed")
Figure 4: Show Proportion Zero Cutoff on Distributions
The chosen cutoff excludes the peak Proportion Zeros from all groups (Figure 4).
Use the cutoff to filter the remaining genes. Only 4k or ~11% of genes pass with this cutoff. Filter the SCE object to this set of genes.
## get a list of filtered genes
filtered_genes <- filter_prop_zero(prop_zeros, cutoff = propZero_limit)
## How many genes pass the filter?
length(filtered_genes)
#> [1] 4005
## What % of genes is this
length(filtered_genes) / nrow(sce.dlpfc.tran)
#> [1] 0.2388336
## Filter the sce object
sce.dlpfc.tran <- sce.dlpfc.tran[filtered_genes, ]One last check of the sparsity, more non-blue means more non-zero values for the Rank Invariance calculation, which prevents rank ties (Figure 5).
## this data down to 50% sparse
sum(assays(sce.dlpfc.tran)$counts == 0) / (nrow(sce.dlpfc.tran) * ncol(sce.dlpfc.tran))
#> [1] 0.5010587
## re-plot heatmap
count_test <- as.matrix(assays(sce.dlpfc.tran)$logcounts[seq_len(1000), seq_len(500)])
pheatmap::pheatmap(count_test,
    cluster_rows = FALSE,
    cluster_cols = FALSE,
    show_rownames = FALSE,
    show_colnames = FALSE
)
Figure 5: Heatmap of the snRNA-seq counts after Proption Zero filtering
To get the Rank Invariance (RI), the rank of the genes across the cells in a group, and between groups is considered. One way to calculate RI is to find the group rank values, and cell rank values separately, then combine them as shown below. The genes with the top RI values are the best candidate TREGs.
## Get the rank of the gene in each group
group_rank <- rank_group(sce.dlpfc.tran, group_col = "cellType.broad")
## Get the rank of the gene for each cell
cell_rank <- rank_cells(sce.dlpfc.tran, group_col = "cellType.broad")
## Use both rankings to calculate rank_invariance()
rank_invar <- rank_invariance(group_rank, cell_rank)
## The top 5 Canidate TREGs:
head(sort(rank_invar, decreasing = TRUE))
#>  MALAT1  CTNND2     FTX    RERE FAM155A    GNAQ 
#>    4005    4004    4003    4002    4001    4000The rank_invariance_express() function combines these three steps into one
function, and achieves the same results.
## rank_invariance_express() runs the previous functions for you
rank_invar2 <- rank_invariance_express(
    sce.dlpfc.tran,
    group_col = "cellType.broad"
)
## Again the top 5 Canidate TREGs:
head(sort(rank_invar2, decreasing = TRUE))
#>  MALAT1  CTNND2     FTX    RERE FAM155A    GNAQ 
#>    4005    4004    4003    4002    4001    4000
## Check computationally that the results are identical
stopifnot(identical(rank_invar, rank_invar2))The rank_invariance_express() function is more efficient as well, as it loops
through the data to rank the genes over cells, and groups at the same time.
We have identified top candidate TREG genes from this dataset, by applying
Proportion Zero filtering and calculating the Rank Invariance using the TREG package.
This provides a list of candidate genes that can be useful for estimating total
RNA expression in assays such as smFISH.
However, we are unable to assess other important qualities of these genes that ensure they are experimentally compatible with the chosen assay. For example, in smFISH with RNAscope it is important that a TREG be expressed at a level that individual puncta can be accurately counted, and have a dynamic range of puncta. During experimental validation we found that MALAT1 was too highly expressed in the human DLPFC to segment individual puncta, and ruled it out as an experimentally useful TREG (Huuki-Myers, Montgomery, Kwon et al., 2022).
Therefore, we recommend that TREGs be evaluated in the assay or analysis of choice you perform a validation experiment with a pilot sample before implementing experiments using it on rare and valuable samples.
If you are designing a sc/snRNA-seq study to use as a reference for deconvolution of bulk RNA-seq, we recommend that you generate spatially-adjacent dissections in order to use them for RNAscope experiments. By doing so, you could identify cell types in your sc/snRNA-seq data, then identify candidate TREGs based on those cell types, and use these candidate TREGs in your spatially-adjacent dissections to quantify size and total RNA amounts for the cell types of interest (Huuki-Myers, Montgomery, Kwon et al., 2022). Furthermore, the RNAscope data alone can be used as a gold standard reference for cell fractions.
TREGs could be useful for other research purposes and other contexts than the ones we envisioned!
Thanks for your interest in TREGs :)
The TREG package (Huuki-Myers and Collado-Torres, 2022) was made possible thanks to:
Code for creating the vignette
## Create the vignette
library("rmarkdown")
system.time(render("finding_Total_RNA_Expression_Genes.Rmd"))
## Extract the R code
library("knitr")
knit("finding_Total_RNA_Expression_Genes.Rmd", tangle = TRUE)Date the vignette was generated.
#> [1] "2022-11-01 19:16:08 EDT"Wallclock time spent generating the vignette.
#> Time difference of 1.312 minsR session information.
#> ─ Session info ───────────────────────────────────────────────────────────────────────────────────────────────────────
#>  setting  value
#>  version  R version 4.2.1 (2022-06-23)
#>  os       Ubuntu 20.04.5 LTS
#>  system   x86_64, linux-gnu
#>  ui       X11
#>  language (EN)
#>  collate  C
#>  ctype    en_US.UTF-8
#>  tz       America/New_York
#>  date     2022-11-01
#>  pandoc   2.5 @ /usr/bin/ (via rmarkdown)
#> 
#> ─ Packages ───────────────────────────────────────────────────────────────────────────────────────────────────────────
#>  package              * version  date (UTC) lib source
#>  assertthat             0.2.1    2019-03-21 [2] CRAN (R 4.2.1)
#>  backports              1.4.1    2021-12-13 [2] CRAN (R 4.2.1)
#>  bibtex                 0.5.0    2022-09-25 [2] CRAN (R 4.2.1)
#>  Biobase              * 2.58.0   2022-11-01 [2] Bioconductor
#>  BiocFileCache          2.6.0    2022-11-01 [2] Bioconductor
#>  BiocGenerics         * 0.44.0   2022-11-01 [2] Bioconductor
#>  BiocManager            1.30.19  2022-10-25 [2] CRAN (R 4.2.1)
#>  BiocStyle            * 2.26.0   2022-11-01 [2] Bioconductor
#>  bit                    4.0.4    2020-08-04 [2] CRAN (R 4.2.1)
#>  bit64                  4.0.5    2020-08-30 [2] CRAN (R 4.2.1)
#>  bitops                 1.0-7    2021-04-24 [2] CRAN (R 4.2.1)
#>  blob                   1.2.3    2022-04-10 [2] CRAN (R 4.2.1)
#>  bookdown               0.29     2022-09-12 [2] CRAN (R 4.2.1)
#>  bslib                  0.4.0    2022-07-16 [2] CRAN (R 4.2.1)
#>  cachem                 1.0.6    2021-08-19 [2] CRAN (R 4.2.1)
#>  cli                    3.4.1    2022-09-23 [2] CRAN (R 4.2.1)
#>  colorspace             2.0-3    2022-02-21 [2] CRAN (R 4.2.1)
#>  curl                   4.3.3    2022-10-06 [2] CRAN (R 4.2.1)
#>  DBI                    1.1.3    2022-06-18 [2] CRAN (R 4.2.1)
#>  dbplyr                 2.2.1    2022-06-27 [2] CRAN (R 4.2.1)
#>  DelayedArray           0.24.0   2022-11-01 [2] Bioconductor
#>  digest                 0.6.30   2022-10-18 [2] CRAN (R 4.2.1)
#>  dplyr                * 1.0.10   2022-09-01 [2] CRAN (R 4.2.1)
#>  ellipsis               0.3.2    2021-04-29 [2] CRAN (R 4.2.1)
#>  evaluate               0.17     2022-10-07 [2] CRAN (R 4.2.1)
#>  fansi                  1.0.3    2022-03-24 [2] CRAN (R 4.2.1)
#>  farver                 2.1.1    2022-07-06 [2] CRAN (R 4.2.1)
#>  fastmap                1.1.0    2021-01-25 [2] CRAN (R 4.2.1)
#>  filelock               1.0.2    2018-10-05 [2] CRAN (R 4.2.1)
#>  generics               0.1.3    2022-07-05 [2] CRAN (R 4.2.1)
#>  GenomeInfoDb         * 1.34.0   2022-11-01 [2] Bioconductor
#>  GenomeInfoDbData       1.2.9    2022-09-28 [2] Bioconductor
#>  GenomicRanges        * 1.50.0   2022-11-01 [2] Bioconductor
#>  ggplot2              * 3.3.6    2022-05-03 [2] CRAN (R 4.2.1)
#>  glue                   1.6.2    2022-02-24 [2] CRAN (R 4.2.1)
#>  gtable                 0.3.1    2022-09-01 [2] CRAN (R 4.2.1)
#>  highr                  0.9      2021-04-16 [2] CRAN (R 4.2.1)
#>  htmltools              0.5.3    2022-07-18 [2] CRAN (R 4.2.1)
#>  httr                   1.4.4    2022-08-17 [2] CRAN (R 4.2.1)
#>  IRanges              * 2.32.0   2022-11-01 [2] Bioconductor
#>  jquerylib              0.1.4    2021-04-26 [2] CRAN (R 4.2.1)
#>  jsonlite               1.8.3    2022-10-21 [2] CRAN (R 4.2.1)
#>  knitr                  1.40     2022-08-24 [2] CRAN (R 4.2.1)
#>  labeling               0.4.2    2020-10-20 [2] CRAN (R 4.2.1)
#>  lattice                0.20-45  2021-09-22 [2] CRAN (R 4.2.1)
#>  lifecycle              1.0.3    2022-10-07 [2] CRAN (R 4.2.1)
#>  lubridate              1.8.0    2021-10-07 [2] CRAN (R 4.2.1)
#>  magick                 2.7.3    2021-08-18 [2] CRAN (R 4.2.1)
#>  magrittr               2.0.3    2022-03-30 [2] CRAN (R 4.2.1)
#>  Matrix                 1.5-1    2022-09-13 [2] CRAN (R 4.2.1)
#>  MatrixGenerics       * 1.10.0   2022-11-01 [2] Bioconductor
#>  matrixStats          * 0.62.0   2022-04-19 [2] CRAN (R 4.2.1)
#>  memoise                2.0.1    2021-11-26 [2] CRAN (R 4.2.1)
#>  munsell                0.5.0    2018-06-12 [2] CRAN (R 4.2.1)
#>  pheatmap             * 1.0.12   2019-01-04 [2] CRAN (R 4.2.1)
#>  pillar                 1.8.1    2022-08-19 [2] CRAN (R 4.2.1)
#>  pkgconfig              2.0.3    2019-09-22 [2] CRAN (R 4.2.1)
#>  plyr                   1.8.7    2022-03-24 [2] CRAN (R 4.2.1)
#>  purrr                  0.3.5    2022-10-06 [2] CRAN (R 4.2.1)
#>  R6                     2.5.1    2021-08-19 [2] CRAN (R 4.2.1)
#>  rafalib                1.0.0    2015-08-09 [2] CRAN (R 4.2.1)
#>  rappdirs               0.3.3    2021-01-31 [2] CRAN (R 4.2.1)
#>  RColorBrewer           1.1-3    2022-04-03 [2] CRAN (R 4.2.1)
#>  Rcpp                   1.0.9    2022-07-08 [2] CRAN (R 4.2.1)
#>  RCurl                  1.98-1.9 2022-10-03 [2] CRAN (R 4.2.1)
#>  RefManageR           * 1.4.0    2022-09-30 [2] CRAN (R 4.2.1)
#>  rlang                  1.0.6    2022-09-24 [2] CRAN (R 4.2.1)
#>  rmarkdown              2.17     2022-10-07 [2] CRAN (R 4.2.1)
#>  RSQLite                2.2.18   2022-10-04 [2] CRAN (R 4.2.1)
#>  S4Vectors            * 0.36.0   2022-11-01 [2] Bioconductor
#>  sass                   0.4.2    2022-07-16 [2] CRAN (R 4.2.1)
#>  scales                 1.2.1    2022-08-20 [2] CRAN (R 4.2.1)
#>  sessioninfo          * 1.2.2    2021-12-06 [2] CRAN (R 4.2.1)
#>  SingleCellExperiment * 1.20.0   2022-11-01 [2] Bioconductor
#>  stringi                1.7.8    2022-07-11 [2] CRAN (R 4.2.1)
#>  stringr                1.4.1    2022-08-20 [2] CRAN (R 4.2.1)
#>  SummarizedExperiment * 1.28.0   2022-11-01 [2] Bioconductor
#>  tibble               * 3.1.8    2022-07-22 [2] CRAN (R 4.2.1)
#>  tidyr                * 1.2.1    2022-09-08 [2] CRAN (R 4.2.1)
#>  tidyselect             1.2.0    2022-10-10 [2] CRAN (R 4.2.1)
#>  TREG                 * 1.2.0    2022-11-01 [1] Bioconductor
#>  utf8                   1.2.2    2021-07-24 [2] CRAN (R 4.2.1)
#>  vctrs                  0.5.0    2022-10-22 [2] CRAN (R 4.2.1)
#>  withr                  2.5.0    2022-03-03 [2] CRAN (R 4.2.1)
#>  xfun                   0.34     2022-10-18 [2] CRAN (R 4.2.1)
#>  xml2                   1.3.3    2021-11-30 [2] CRAN (R 4.2.1)
#>  XVector                0.38.0   2022-11-01 [2] Bioconductor
#>  yaml                   2.3.6    2022-10-18 [2] CRAN (R 4.2.1)
#>  zlibbioc               1.44.0   2022-11-01 [2] Bioconductor
#> 
#>  [1] /tmp/RtmpdXOyw6/Rinst27cff4bc6c1e1
#>  [2] /home/biocbuild/bbs-3.16-bioc/R/library
#> 
#> ──────────────────────────────────────────────────────────────────────────────────────────────────────────────────────This vignette was generated using BiocStyle (Oleś, 2022), knitr (Xie, 2022) and rmarkdown (Allaire, Xie, McPherson et al., 2022) running behind the scenes.
Citations made with RefManageR (McLean, 2017).
[1] J. Allaire, Y. Xie, J. McPherson, et al. rmarkdown: Dynamic Documents for R. R package version 2.17. 2022. URL: https://github.com/rstudio/rmarkdown.
[2] D. Bates, M. Maechler, and M. Jagan. Matrix: Sparse and Dense Matrix Classes and Methods. R package version 1.5-1. 2022. URL: https://CRAN.R-project.org/package=Matrix.
[3] L. Henry and H. Wickham. purrr: Functional Programming Tools. R package version 0.3.5. 2022. URL: https://CRAN.R-project.org/package=purrr.
[4] L. A. Huuki-Myers and L. Collado-Torres. TREG: a R/Bioconductor package to identify Total RNA Expression Genes. https://github.com/LieberInstitute/TREG/TREG - R package version 1.2.0. 2022. DOI: 10.18129/B9.bioc.TREG. URL: http://www.bioconductor.org/packages/TREG.
[5] L. A. Huuki-Myers, K. D. Montgomery, S. H. Kwon, et al. “Data Driven Identification of Total RNA Expression Genes”TREGs" for estimation of RNA abundance in heterogeneous cell types”. In: bioRxiv (2022). DOI: 10.1101/2022.04.28.489923. URL: https://doi.org/10.1101/2022.04.28.489923.
[6] R. A. Irizarry and M. I. Love. rafalib: Convenience Functions for Routine Data Exploration. R package version 1.0.0. 2015. URL: https://CRAN.R-project.org/package=rafalib.
[7] R. Kolde. pheatmap: Pretty Heatmaps. R package version 1.0.12. 2019. URL: https://CRAN.R-project.org/package=pheatmap.
[8] 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.
[9] M. Morgan, V. Obenchain, J. Hester, et al. SummarizedExperiment: SummarizedExperiment container. R package version 1.28.0. 2022. URL: https://bioconductor.org/packages/SummarizedExperiment.
[10] K. Müller and H. Wickham. tibble: Simple Data Frames. R package version 3.1.8. 2022. URL: https://CRAN.R-project.org/package=tibble.
[11] A. Oleś. BiocStyle: Standard styles for vignettes and other Bioconductor documents. R package version 2.26.0. 2022. URL: https://github.com/Bioconductor/BiocStyle.
[12] R Core Team. R: A Language and Environment for Statistical Computing. R Foundation for Statistical Computing. Vienna, Austria, 2022. URL: https://www.R-project.org/.
[13] L. Shepherd and M. Morgan. BiocFileCache: Manage Files Across Sessions. R package version 2.6.0. 2022.
[14] M. N. Tran, K. R. Maynard, A. Spangler, et al. “Single-nucleus transcriptome analysis reveals cell-type-specific molecular signatures across reward circuitry in the human brain”. In: Neuron (2021). DOI: 10.1016/j.neuron.2021.09.001.
[15] H. Wickham. ggplot2: Elegant Graphics for Data Analysis. Springer-Verlag New York, 2016. ISBN: 978-3-319-24277-4. URL: https://ggplot2.tidyverse.org.
[16] 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.
[17] H. Wickham, W. Chang, R. Flight, et al. sessioninfo: R Session Information. R package version 1.2.2. 2021. URL: https://CRAN.R-project.org/package=sessioninfo.
[18] H. Wickham, R. François, L. Henry, et al. dplyr: A Grammar of Data Manipulation. R package version 1.0.10. 2022. URL: https://CRAN.R-project.org/package=dplyr.
[19] H. Wickham and M. Girlich. tidyr: Tidy Messy Data. R package version 1.2.1. 2022. URL: https://CRAN.R-project.org/package=tidyr.
[20] Y. Xie. knitr: A General-Purpose Package for Dynamic Report Generation in R. R package version 1.40. 2022. URL: https://yihui.org/knitr/.