Contents

1 Introduction

nnSVG is a method for scalable identification of spatially variable genes (SVGs) in spatially-resolved transcriptomics data.

The nnSVG method is based on nearest-neighbor Gaussian processes (Datta et al., 2016, Finley et al., 2019) and uses the BRISC algorithm (Saha and Datta, 2018) for model fitting and parameter estimation. nnSVG allows identification and ranking of SVGs with flexible length scales across a tissue slide or within spatial domains defined by covariates. The method scales linearly with the number of spatial locations and can be applied to datasets containing thousands or more spatial locations.

nnSVG is implemented as an R package within the Bioconductor framework, and is available from Bioconductor.

More details describing the method are available in our preprint, available from bioRxiv.

2 Installation

The following code will install the latest release version of the nnSVG package from Bioconductor. Additional details are shown on the Bioconductor page.

if (!requireNamespace("BiocManager", quietly = TRUE)) {
  install.packages("BiocManager")
}
BiocManager::install("nnSVG")

The latest development version can also be installed from the devel version of Bioconductor or from GitHub.

3 Input data format

In the examples below, we assume the input data are provided as a SpatialExperiment Bioconductor object. In this case, the outputs are stored in the rowData of the SpatialExperiment object.

However, the inputs can also be provided as a numeric matrix of normalized and transformed counts (e.g. log-transformed normalized counts, also known as logcounts) and a numeric matrix of spatial coordinates.

To provide the inputs as numeric matrices, please install the development version of the package from GitHub or the devel version of Bioconductor (which will become the new Bioconductor release version in October 2022).

4 Tutorial

4.1 Standard analysis

4.1.1 Run nnSVG

Here we show a short example demonstrating how to run nnSVG.

For faster runtime in this example, we subsample the dataset and run nnSVG on only a small number of genes. For a full analysis, the subsampling step can be skipped.

library(SpatialExperiment)
library(STexampleData)
library(scran)
library(nnSVG)
library(ggplot2)
# load example dataset from STexampleData package
spe <- Visium_humanDLPFC()
## see ?STexampleData and browseVignettes('STexampleData') for documentation
## loading from cache
dim(spe)
## [1] 33538  4992
# preprocessing steps

# keep only spots over tissue
spe <- spe[, colData(spe)$in_tissue == 1]
dim(spe)
## [1] 33538  3639
# skip spot-level quality control, since this has been performed previously 
# on this dataset

# filter low-expressed and mitochondrial genes
# using default filtering parameters
spe <- filter_genes(spe)
## Gene filtering: removing mitochondrial genes
## removed 13 mitochondrial genes
## Gene filtering: retaining genes with at least 3 counts in at least 0.5% (n = 19) of spatial locations
## removed 30216 out of 33525 genes due to low expression
# calculate logcounts (log-transformed normalized counts) using scran package
# using library size factors
spe <- computeLibraryFactors(spe)
spe <- logNormCounts(spe)

assayNames(spe)
## [1] "counts"    "logcounts"
# select small set of random genes and several known SVGs for 
# faster runtime in this example
set.seed(123)
ix_random <- sample(seq_len(nrow(spe)), 10)

known_genes <- c("MOBP", "PCP4", "SNAP25", "HBB", "IGKC", "NPY")
ix_known <- which(rowData(spe)$gene_name %in% known_genes)

ix <- c(ix_known, ix_random)

spe <- spe[ix, ]
dim(spe)
## [1]   16 3639
# run nnSVG

# set seed for reproducibility
set.seed(123)
# using a single thread in this example
spe <- nnSVG(spe)

# show results
rowData(spe)
## DataFrame with 16 rows and 17 columns
##                         gene_id   gene_name    feature_type   sigma.sq
##                     <character> <character>     <character>  <numeric>
## ENSG00000211592 ENSG00000211592        IGKC Gene Expression   0.565654
## ENSG00000168314 ENSG00000168314        MOBP Gene Expression   1.387394
## ENSG00000122585 ENSG00000122585         NPY Gene Expression   0.285674
## ENSG00000244734 ENSG00000244734         HBB Gene Expression   0.329421
## ENSG00000132639 ENSG00000132639      SNAP25 Gene Expression   0.430040
## ...                         ...         ...             ...        ...
## ENSG00000130382 ENSG00000130382       MLLT1 Gene Expression 0.00978555
## ENSG00000036672 ENSG00000036672        USP2 Gene Expression 0.00307277
## ENSG00000086232 ENSG00000086232     EIF2AK1 Gene Expression 0.00315782
## ENSG00000106278 ENSG00000106278      PTPRZ1 Gene Expression 0.00279851
## ENSG00000133606 ENSG00000133606       MKRN1 Gene Expression 0.00632245
##                    tau.sq        phi    loglik   runtime      mean       var
##                 <numeric>  <numeric> <numeric> <numeric> <numeric> <numeric>
## ENSG00000211592  0.455041   20.10688  -4531.64     1.651  0.622937  1.007454
## ENSG00000168314  0.364188    1.10202  -3663.60     0.937  0.805525  1.205673
## ENSG00000122585  0.280173   71.65329  -3995.23     1.352  0.393975  0.567383
## ENSG00000244734  0.353754   27.81410  -4044.96     2.332  0.411262  0.697673
## ENSG00000132639  0.430106    3.03385  -3912.70     1.013  3.451926  0.857922
## ...                   ...        ...       ...       ...       ...       ...
## ENSG00000130382  0.283115 50.9880748  -2927.61     1.099  0.298698  0.292976
## ENSG00000036672  0.241105 12.5382833  -2597.00     0.793  0.248384  0.244218
## ENSG00000086232  0.266973 25.9302215  -2781.47     1.384  0.275193  0.270208
## ENSG00000106278  0.367893  9.5280046  -3357.32     1.584  0.352159  0.370784
## ENSG00000133606  0.272432  0.0827087  -2831.51     1.874  0.295404  0.278806
##                     spcov    prop_sv loglik_lm   LR_stat      rank        pval
##                 <numeric>  <numeric> <numeric> <numeric> <numeric>   <numeric>
## ENSG00000211592  1.207346   0.554185  -5176.53  1289.775         3           0
## ENSG00000168314  1.462248   0.792080  -5503.33  3679.464         1           0
## ENSG00000122585  1.356646   0.504861  -4131.87   273.278         6           0
## ENSG00000244734  1.395587   0.482191  -4507.99   926.046         4           0
## ENSG00000132639  0.189973   0.499961  -4884.19  1942.986         2           0
## ...                   ...        ...       ...       ...       ...         ...
## ENSG00000130382  0.331177 0.03340915  -2929.28   3.35216        12 0.187106089
## ENSG00000036672  0.223173 0.01258414  -2598.08   2.15483        13 0.340473854
## ENSG00000086232  0.204200 0.01168999  -2782.09   1.23716        14 0.538708116
## ENSG00000106278  0.150219 0.00754943  -3357.83   1.01111        15 0.603169436
## ENSG00000133606  0.269170 0.02268111  -2839.08  15.15227         9 0.000512539
##                       padj
##                  <numeric>
## ENSG00000211592          0
## ENSG00000168314          0
## ENSG00000122585          0
## ENSG00000244734          0
## ENSG00000132639          0
## ...                    ...
## ENSG00000130382 0.24947479
## ENSG00000036672 0.41904474
## ENSG00000086232 0.61566642
## ENSG00000106278 0.64338073
## ENSG00000133606 0.00091118

4.1.2 Investigate results

The results are stored in the rowData of the SpatialExperiment object.

The main results of interest are:

  • LR_stat: likelihood ratio (LR) statistics
  • rank: rank of top SVGs according to LR statistics
  • pval: p-values from asymptotic chi-squared distribution with 2 degrees of freedom
  • padj: p-values adjusted for multiple testing, which can be used to define a cutoff for statistically significant SVGs (e.g. padj <= 0.05)
  • prop_sv: effect size, defined as proportion of spatial variance out of total variance
# number of significant SVGs
table(rowData(spe)$padj <= 0.05)
## 
## FALSE  TRUE 
##     7     9
# show results for top n SVGs
rowData(spe)[order(rowData(spe)$rank)[1:10], ]
## DataFrame with 10 rows and 17 columns
##                         gene_id   gene_name    feature_type   sigma.sq
##                     <character> <character>     <character>  <numeric>
## ENSG00000168314 ENSG00000168314        MOBP Gene Expression 1.38739399
## ENSG00000132639 ENSG00000132639      SNAP25 Gene Expression 0.43003959
## ENSG00000211592 ENSG00000211592        IGKC Gene Expression 0.56565436
## ENSG00000244734 ENSG00000244734         HBB Gene Expression 0.32942113
## ENSG00000183036 ENSG00000183036        PCP4 Gene Expression 0.23102221
## ENSG00000122585 ENSG00000122585         NPY Gene Expression 0.28567358
## ENSG00000129562 ENSG00000129562        DAD1 Gene Expression 0.02389606
## ENSG00000114923 ENSG00000114923      SLC4A3 Gene Expression 0.01147170
## ENSG00000133606 ENSG00000133606       MKRN1 Gene Expression 0.00632245
## ENSG00000143543 ENSG00000143543         JTB Gene Expression 0.07547797
##                    tau.sq         phi    loglik   runtime      mean       var
##                 <numeric>   <numeric> <numeric> <numeric> <numeric> <numeric>
## ENSG00000168314  0.364188   1.1020177  -3663.60     0.937  0.805525  1.205673
## ENSG00000132639  0.430106   3.0338473  -3912.70     1.013  3.451926  0.857922
## ENSG00000211592  0.455041  20.1068839  -4531.64     1.651  0.622937  1.007454
## ENSG00000244734  0.353754  27.8140976  -4044.96     2.332  0.411262  0.697673
## ENSG00000183036  0.452735   8.2722785  -4026.22     1.274  0.687961  0.684598
## ENSG00000122585  0.280173  71.6532892  -3995.23     1.352  0.393975  0.567383
## ENSG00000129562  0.464723  10.1418819  -3842.24     1.706  0.549318  0.489167
## ENSG00000114923  0.237260  12.7656826  -2617.36     1.193  0.250768  0.248816
## ENSG00000133606  0.272432   0.0827087  -2831.51     1.874  0.295404  0.278806
## ENSG00000143543  0.463561 119.7470905  -4036.28     1.592  0.654919  0.539172
##                     spcov   prop_sv loglik_lm    LR_stat      rank        pval
##                 <numeric> <numeric> <numeric>  <numeric> <numeric>   <numeric>
## ENSG00000168314  1.462248 0.7920804  -5503.33 3679.46397         1 0.00000e+00
## ENSG00000132639  0.189973 0.4999614  -4884.19 1942.98556         2 0.00000e+00
## ENSG00000211592  1.207346 0.5541853  -5176.53 1289.77508         3 0.00000e+00
## ENSG00000244734  1.395587 0.4821910  -4507.99  926.04573         4 0.00000e+00
## ENSG00000183036  0.698656 0.3378716  -4473.57  894.68884         5 0.00000e+00
## ENSG00000122585  1.356646 0.5048608  -4131.87  273.27818         6 0.00000e+00
## ENSG00000129562  0.281410 0.0489053  -3861.98   39.49098         7 2.65854e-09
## ENSG00000114923  0.427112 0.0461207  -2632.02   29.31376         8 4.31119e-07
## ENSG00000133606  0.269170 0.0226811  -2839.08   15.15227         9 5.12539e-04
## ENSG00000143543  0.419491 0.1400231  -4039.07    5.59669        10 6.09108e-02
##                        padj
##                   <numeric>
## ENSG00000168314 0.00000e+00
## ENSG00000132639 0.00000e+00
## ENSG00000211592 0.00000e+00
## ENSG00000244734 0.00000e+00
## ENSG00000183036 0.00000e+00
## ENSG00000122585 0.00000e+00
## ENSG00000129562 6.07667e-09
## ENSG00000114923 8.62238e-07
## ENSG00000133606 9.11180e-04
## ENSG00000143543 9.74572e-02
# plot spatial expression of top-ranked SVG
ix <- which(rowData(spe)$rank == 1)
ix_name <- rowData(spe)$gene_name[ix]
ix_name
## [1] "MOBP"
df <- as.data.frame(
  cbind(spatialCoords(spe), 
        expr = counts(spe)[ix, ]))

ggplot(df, aes(x = pxl_col_in_fullres, y = pxl_row_in_fullres, color = expr)) + 
  geom_point(size = 0.8) + 
  coord_fixed() + 
  scale_y_reverse() + 
  scale_color_gradient(low = "gray90", high = "blue", 
                       trans = "sqrt", breaks = range(df$expr), 
                       name = "counts") + 
  ggtitle(ix_name) + 
  theme_bw() + 
  theme(plot.title = element_text(face = "italic"), 
        panel.grid = element_blank(), 
        axis.title = element_blank(), 
        axis.text = element_blank(), 
        axis.ticks = element_blank())

4.2 With covariates

4.2.1 Run nnSVG

nnSVG can also be run on datasets with known covariates, such as spatial domains or cell types. In this case, the method will identify SVGs within regions defined by the selected covariates, e.g. within spatial domains.

Here we run a short example demonstrating this functionality.

As above, for faster runtime in this example, we subsample a very small subset of the data. For a full analysis, the subsampling step can be skipped.

library(SpatialExperiment)
library(STexampleData)
library(scran)
library(nnSVG)
library(ggplot2)
# load example dataset from STexampleData package
spe <- SlideSeqV2_mouseHPC()
## see ?STexampleData and browseVignettes('STexampleData') for documentation
## loading from cache
dim(spe)
## [1] 23264 53208
# preprocessing steps

# remove spots with NA cell type labels
spe <- spe[, !is.na(colData(spe)$celltype)]
dim(spe)
## [1] 23264 15003
# check cell type labels
table(colData(spe)$celltype)
## 
##             Astrocyte                   CA1                   CA3 
##                  6688                  1320                  1729 
##         Cajal_Retzius               Choroid                Denate 
##                    60                    21                  1713 
##     Endothelial_Stalk       Endothelial_Tip           Entorihinal 
##                    96                   118                   111 
##             Ependymal           Interneuron Microglia_Macrophages 
##                    94                   730                   299 
##                 Mural          Neurogenesis        Neuron.Slc17a6 
##                    26                    47                    44 
##       Oligodendrocyte        Polydendrocyte 
##                  1793                   114
# filter low-expressed and mitochondrial genes
# using adjusted filtering parameters for this platform
spe <- filter_genes(
  spe, 
  filter_genes_ncounts = 1, 
  filter_genes_pcspots = 1, 
  filter_mito = TRUE
)
## Gene filtering: removing mitochondrial genes
## removed 25 mitochondrial genes
## Gene filtering: retaining genes with at least 1 counts in at least 1% (n = 151) of spatial locations
## removed 14356 out of 23239 genes due to low expression
dim(spe)
## [1]  8883 15003
# calculate log-transformed normalized counts using scran package
# using library size normalization
spe <- computeLibraryFactors(spe)
spe <- logNormCounts(spe)

assayNames(spe)
## [1] "counts"    "logcounts"
# select small set of random genes and several known SVGs for 
# faster runtime in this example
set.seed(123)
ix_random <- sample(seq_len(nrow(spe)), 10)

known_genes <- c("Cpne9", "Rgs14")
ix_known <- which(rowData(spe)$gene_name %in% known_genes)

ix <- c(ix_known, ix_random)

spe <- spe[ix, ]
dim(spe)
## [1]    12 15003
# run nnSVG with covariates

# create model matrix for cell type labels
X <- model.matrix(~ colData(spe)$celltype)
dim(X)
## [1] 15003    17
stopifnot(nrow(X) == ncol(spe))

# set seed for reproducibility
set.seed(123)
# using a single thread in this example
spe <- nnSVG(spe, X = X)

# show results
rowData(spe)
## DataFrame with 12 rows and 15 columns
##          gene_name    sigma.sq    tau.sq         phi    loglik   runtime
##        <character>   <numeric> <numeric>   <numeric> <numeric> <numeric>
## Cpne9        Cpne9 6.91469e-04 0.0275984 1.28476e+01   5484.78    21.985
## Rgs14        Rgs14 2.77038e-03 0.0452680 1.90815e-06   1650.39    16.814
## Eogt          Eogt 8.52208e-05 0.0198481 7.81846e+01   8083.08     9.699
## Ergic1      Ergic1 2.54389e-04 0.0465677 9.77238e+01   1677.16     6.061
## Zfp330      Zfp330 1.21717e-04 0.0266163 1.54158e+02   5879.76    13.103
## ...            ...         ...       ...         ...       ...       ...
## Hdac11      Hdac11 3.54443e-03 0.1648232    145.2971  -7920.47    13.423
## Nab2          Nab2 3.99615e-05 0.0145549     67.5465  10421.32    12.337
## Senp3        Senp3 7.97660e-05 0.0253771     36.8022   6248.29     9.511
## Fbxo6        Fbxo6 9.53698e-05 0.0253767    168.3817   6243.54    11.931
## Ntng1        Ntng1 2.39135e-03 0.0433036     24.1539   1951.83     9.438
##             mean       var     spcov    prop_sv loglik_lm    LR_stat      rank
##        <numeric> <numeric> <numeric>  <numeric> <numeric>  <numeric> <numeric>
## Cpne9  0.0194924 0.0283971  1.349026 0.02444229   5455.66   58.23907         3
## Rgs14  0.0360053 0.0485672  1.461850 0.05767000   1515.16  270.46145         1
## Eogt   0.0125243 0.0200072  0.737089 0.00427529   8082.32    1.52763         6
## Ergic1 0.0359343 0.0471421  0.443854 0.00543309   1676.61    1.10864         8
## Zfp330 0.0201203 0.0267799  0.548330 0.00455223   5879.05    1.42306         7
## ...          ...       ...       ...        ...       ...        ...       ...
## Hdac11 0.1258726 0.1705333  0.472979 0.02105171  -7923.51   6.072553         4
## Nab2   0.0120832 0.0146742  0.523164 0.00273805  10421.82  -1.016745        12
## Senp3  0.0199706 0.0255396  0.447217 0.00313338   6248.60  -0.631638        11
## Fbxo6  0.0188249 0.0255584  0.518768 0.00374409   6243.59  -0.103463        10
## Ntng1  0.0341843 0.0466474  1.430524 0.05233303   1854.40 194.856131         2
##               pval        padj
##          <numeric>   <numeric>
## Cpne9  2.25708e-13 9.02833e-13
## Rgs14  0.00000e+00 0.00000e+00
## Eogt   4.65886e-01 8.41531e-01
## Ergic1 5.74464e-01 8.61696e-01
## Zfp330 4.90893e-01 8.41531e-01
## ...            ...         ...
## Hdac11   0.0480133     0.14404
## Nab2     1.0000000     1.00000
## Senp3    1.0000000     1.00000
## Fbxo6    1.0000000     1.00000
## Ntng1    0.0000000     0.00000

Note that if there are any empty levels in the factor used to create the design matrix, these can be removed as follows:

# create model matrix after dropping empty factor levels
X <- model.matrix(~ droplevels(as.factor(colData(spe)$celltype)))

4.2.2 Investigate results

# number of significant SVGs
table(rowData(spe)$padj <= 0.05)
## 
## FALSE  TRUE 
##     9     3
# show results for top n SVGs
rowData(spe)[order(rowData(spe)$rank)[1:10], ]
## DataFrame with 10 rows and 15 columns
##          gene_name    sigma.sq    tau.sq         phi    loglik   runtime
##        <character>   <numeric> <numeric>   <numeric> <numeric> <numeric>
## Rgs14        Rgs14 2.77038e-03 0.0452680 1.90815e-06   1650.39    16.814
## Ntng1        Ntng1 2.39135e-03 0.0433036 2.41539e+01   1951.83     9.438
## Cpne9        Cpne9 6.91469e-04 0.0275984 1.28476e+01   5484.78    21.985
## Hdac11      Hdac11 3.54443e-03 0.1648232 1.45297e+02  -7920.47    13.423
## Cul1          Cul1 4.19371e-04 0.0688365 9.03178e+01  -1259.18    12.649
## Eogt          Eogt 8.52208e-05 0.0198481 7.81846e+01   8083.08     9.699
## Zfp330      Zfp330 1.21717e-04 0.0266163 1.54158e+02   5879.76    13.103
## Ergic1      Ergic1 2.54389e-04 0.0465677 9.77238e+01   1677.16     6.061
## Gclc          Gclc 1.22209e-04 0.0345001 1.19972e+02   3941.26     7.484
## Fbxo6        Fbxo6 9.53698e-05 0.0253767 1.68382e+02   6243.54    11.931
##             mean       var     spcov    prop_sv loglik_lm     LR_stat      rank
##        <numeric> <numeric> <numeric>  <numeric> <numeric>   <numeric> <numeric>
## Rgs14  0.0360053 0.0485672  1.461850 0.05767000   1515.16 270.4614462         1
## Ntng1  0.0341843 0.0466474  1.430524 0.05233303   1854.40 194.8561307         2
## Cpne9  0.0194924 0.0283971  1.349026 0.02444229   5455.66  58.2390704         3
## Hdac11 0.1258726 0.1705333  0.472979 0.02105171  -7923.51   6.0725529         4
## Cul1   0.0548957 0.0696223  0.373045 0.00605538  -1260.02   1.6933490         5
## Eogt   0.0125243 0.0200072  0.737089 0.00427529   8082.32   1.5276268         6
## Zfp330 0.0201203 0.0267799  0.548330 0.00455223   5879.05   1.4230584         7
## Ergic1 0.0359343 0.0471421  0.443854 0.00543309   1676.61   1.1086351         8
## Gclc   0.0256188 0.0347265  0.431513 0.00352979   3941.30  -0.0749267         9
## Fbxo6  0.0188249 0.0255584  0.518768 0.00374409   6243.59  -0.1034630        10
##               pval        padj
##          <numeric>   <numeric>
## Rgs14  0.00000e+00 0.00000e+00
## Ntng1  0.00000e+00 0.00000e+00
## Cpne9  2.25708e-13 9.02833e-13
## Hdac11 4.80133e-02 1.44040e-01
## Cul1   4.28839e-01 8.41531e-01
## Eogt   4.65886e-01 8.41531e-01
## Zfp330 4.90893e-01 8.41531e-01
## Ergic1 5.74464e-01 8.61696e-01
## Gclc   1.00000e+00 1.00000e+00
## Fbxo6  1.00000e+00 1.00000e+00
# plot spatial expression of top-ranked SVG
ix <- which(rowData(spe)$rank == 1)
ix_name <- rowData(spe)$gene_name[ix]
ix_name
## [1] "Rgs14"
df <- as.data.frame(
  cbind(spatialCoords(spe), 
        expr = counts(spe)[ix, ]))

ggplot(df, aes(x = xcoord, y = ycoord, color = expr)) + 
  geom_point(size = 0.1) + 
  coord_fixed() + 
  scale_color_gradient(low = "gray90", high = "blue", 
                       trans = "sqrt", breaks = range(df$expr), 
                       name = "counts") + 
  ggtitle(ix_name) + 
  theme_bw() + 
  theme(plot.title = element_text(face = "italic"), 
        panel.grid = element_blank(), 
        axis.title = element_blank(), 
        axis.text = element_blank(), 
        axis.ticks = element_blank())

5 Session information

sessionInfo()
## R version 4.2.1 (2022-06-23)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 20.04.5 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.16-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.16-bioc/R/lib/libRlapack.so
## 
## locale:
##  [1] LC_CTYPE=en_US.UTF-8       LC_NUMERIC=C              
##  [3] LC_TIME=en_GB              LC_COLLATE=C              
##  [5] LC_MONETARY=en_US.UTF-8    LC_MESSAGES=en_US.UTF-8   
##  [7] LC_PAPER=en_US.UTF-8       LC_NAME=C                 
##  [9] LC_ADDRESS=C               LC_TELEPHONE=C            
## [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C       
## 
## attached base packages:
## [1] stats4    stats     graphics  grDevices utils     datasets  methods  
## [8] base     
## 
## other attached packages:
##  [1] ggplot2_3.3.6               nnSVG_1.2.0                
##  [3] scran_1.26.0                scuttle_1.8.0              
##  [5] STexampleData_1.5.2         ExperimentHub_2.6.0        
##  [7] AnnotationHub_3.6.0         BiocFileCache_2.6.0        
##  [9] dbplyr_2.2.1                SpatialExperiment_1.8.0    
## [11] SingleCellExperiment_1.20.0 SummarizedExperiment_1.28.0
## [13] Biobase_2.58.0              GenomicRanges_1.50.0       
## [15] GenomeInfoDb_1.34.0         IRanges_2.32.0             
## [17] S4Vectors_0.36.0            BiocGenerics_0.44.0        
## [19] MatrixGenerics_1.10.0       matrixStats_0.62.0         
## [21] BiocStyle_2.26.0           
## 
## loaded via a namespace (and not attached):
##   [1] colorspace_2.0-3              rjson_0.2.21                 
##   [3] ellipsis_0.3.2                bluster_1.8.0                
##   [5] XVector_0.38.0                BiocNeighbors_1.16.0         
##   [7] farver_2.1.1                  bit64_4.0.5                  
##   [9] interactiveDisplayBase_1.36.0 AnnotationDbi_1.60.0         
##  [11] fansi_1.0.3                   codetools_0.2-18             
##  [13] R.methodsS3_1.8.2             sparseMatrixStats_1.10.0     
##  [15] cachem_1.0.6                  knitr_1.40                   
##  [17] jsonlite_1.8.3                cluster_2.1.4                
##  [19] png_0.1-7                     R.oo_1.25.0                  
##  [21] shiny_1.7.3                   HDF5Array_1.26.0             
##  [23] BiocManager_1.30.19           compiler_4.2.1               
##  [25] httr_1.4.4                    dqrng_0.3.0                  
##  [27] assertthat_0.2.1              Matrix_1.5-1                 
##  [29] fastmap_1.1.0                 limma_3.54.0                 
##  [31] cli_3.4.1                     later_1.3.0                  
##  [33] BiocSingular_1.14.0           BRISC_1.0.5                  
##  [35] htmltools_0.5.3               tools_4.2.1                  
##  [37] rsvd_1.0.5                    igraph_1.3.5                 
##  [39] gtable_0.3.1                  glue_1.6.2                   
##  [41] GenomeInfoDbData_1.2.9        RANN_2.6.1                   
##  [43] dplyr_1.0.10                  rappdirs_0.3.3               
##  [45] Rcpp_1.0.9                    jquerylib_0.1.4              
##  [47] vctrs_0.5.0                   Biostrings_2.66.0            
##  [49] rhdf5filters_1.10.0           DelayedMatrixStats_1.20.0    
##  [51] xfun_0.34                     stringr_1.4.1                
##  [53] beachmat_2.14.0               mime_0.12                    
##  [55] lifecycle_1.0.3               irlba_2.3.5.1                
##  [57] statmod_1.4.37                edgeR_3.40.0                 
##  [59] scales_1.2.1                  zlibbioc_1.44.0              
##  [61] promises_1.2.0.1              parallel_4.2.1               
##  [63] rhdf5_2.42.0                  yaml_2.3.6                   
##  [65] curl_4.3.3                    pbapply_1.5-0                
##  [67] memoise_2.0.1                 sass_0.4.2                   
##  [69] stringi_1.7.8                 RSQLite_2.2.18               
##  [71] highr_0.9                     BiocVersion_3.16.0           
##  [73] ScaledMatrix_1.6.0            filelock_1.0.2               
##  [75] rdist_0.0.5                   BiocParallel_1.32.0          
##  [77] rlang_1.0.6                   pkgconfig_2.0.3              
##  [79] bitops_1.0-7                  evaluate_0.17                
##  [81] lattice_0.20-45               purrr_0.3.5                  
##  [83] Rhdf5lib_1.20.0               labeling_0.4.2               
##  [85] bit_4.0.4                     tidyselect_1.2.0             
##  [87] magrittr_2.0.3                bookdown_0.29                
##  [89] R6_2.5.1                      magick_2.7.3                 
##  [91] generics_0.1.3                metapod_1.6.0                
##  [93] DelayedArray_0.24.0           DBI_1.1.3                    
##  [95] pillar_1.8.1                  withr_2.5.0                  
##  [97] KEGGREST_1.38.0               RCurl_1.98-1.9               
##  [99] tibble_3.1.8                  crayon_1.5.2                 
## [101] DropletUtils_1.18.0           utf8_1.2.2                   
## [103] rmarkdown_2.17                locfit_1.5-9.6               
## [105] grid_4.2.1                    blob_1.2.3                   
## [107] digest_0.6.30                 xtable_1.8-4                 
## [109] httpuv_1.6.6                  R.utils_2.12.1               
## [111] munsell_0.5.0                 bslib_0.4.0