The package spatialGE provides a collection of tools for the visualization of gene expression from spatially-resolved transcriptomic experiments. The data input methods have been designed so that any data can be analyzed as long as it contains gene expression counts per region of interest (ROI), spot, or cell, and the spatial coordinates of those ROIs, spots, or cells, as it is generated in platforms such as GeoMx, Visium, and CosMx-SMI. Specialized data input options are available to allow interoperability with workflows such as Space Ranger.

Installation

The spatialGE repository is available at GitHub and can be installed via devtools. To install devtools (in case it is not already installed in your R), please run the following code:

if("devtools" %in% rownames(installed.packages()) == FALSE){
  install.packages("devtools")
}

After making sure devtools is installed, proceed to install spatialGE:

#devtools::install_github("fridleylab/spatialGE")

To use spatialGE, load the package.

library('spatialGE')

Spatially-resolved expression of triple negative breast cancer tumor biopsies

To show the utility of some of the functions in spatialGE, we use the spatial transcriptomics data set generated with the platform Visium by Bassiouni et al. (2023). This data set includes triple negative breast cancer biopsies from 22 patients, with two tissue slices per patient. The Visium platform allows gene expression quantitation in approximately 5000 capture areas (i.e., “spots”) separated each other by 100μM. The spots are 55μM in diameter, which corresponds to 1-10 cells according to the manufacturer.

Data for all the tissue slices are available at the Gene Expression Omnibus (GEO). For the purpose of this tutorial, we will use eight samples from four patients. The GEO repositories can be accessed using the following links:

Within each link, locate the files filtered_feature_bc_matrix.h5, tissue_positions_list.csv.gz, tissue_hires_image.png.gz, and scalefactors_json.json.gz files and download them. Create a folder in your Desktop with the name tnbc_visium. Within that folder, create eight folders, each with the sample names previously shown (i.e., sample_117d, sample_117e, sample_120d, …). Next, place within each folder the corresponding filtered_feature_bc_matrix.h5 file. Then, create within each folder, another folder named spatial and place the corresponding tissue_positions_list.csv.gz, tissue_hires_image.png.gz, and scalefactors_json.json.gz files within that spatial folder. The files within spatial must be decompressed, which is achieved by double-clicking each file. This file structure, corresponds roughly to the file structure in outputs generated by Space Ranger, the software provided by 10X Genomics to process Visium data.

Creating an STList (Spatial Transcriptomics List)

In spatialGE, raw and processed data are stored in an STlist (S4 class object). The STlist can be created with the function STlist, which can take data using different formats (see here for more info or type ?STlist in the R console). In this tutorial we will provide the file paths to the folders created in the previous step.

Additionally, we will input meta data associated with each sample. The meta data is provided in the form of a comma-delimited file. We have extracted some of the clinical meta data for the eight samples from the original publication and saved it as part of the spatialGE package. The user is encouraged to look at the structure of this file by downloading it from the GitHub repository. The most important aspect when constructing this file is that the sample names are in the first column, and they match the names of the folders containing the data:

Once the folders with samples have been created, we can use R to generate the file paths to be passed to the STlist function. Assuming that the user created the folder in the computer’s Desktop:

visium_folders <- list.dirs('~/Desktop/tnbc_visium/', full.names=T, recursive=F)

The meta data can be accessed directly from the spatialGE package installed in the computer like so:

data_files <- system.file("extdata/tnbc_bassiouni", package="spatialGE")
clin_file <- list.files(data_files, full.names=T, recursive=F, pattern='clinical')

We can load the files into an STlist using this command:

tnbc <- STlist(rnacounts=visium_folders, samples=clin_file)
#> Found Visium data
#> Found 8 Visium samples
#>   Data read completed
#> Matching gene expression and coordinate data...
#> Completed STlist!
#> 

The tnbc object is an STlist that contains the count data, spot coordinates, and clinical meta data.

tnbc
#> Spatial Transcriptomics List (STlist).
#> 8 spatial array(s):
#>  sample_117d (1949 ROIs|spots|cells x 23755 genes)
#>  sample_117e (754 ROIs|spots|cells x 21662 genes)
#>  sample_120d (1359 ROIs|spots|cells x 21886 genes)
#>  sample_120e (844 ROIs|spots|cells x 20389 genes)
#>  sample_092a (1289 ROIs|spots|cells x 24227 genes)
#>  sample_092b (1376 ROIs|spots|cells x 23511 genes)
#>  sample_094c (1542 ROIs|spots|cells x 24581 genes)
#>  sample_094d (1391 ROIs|spots|cells x 24167 genes)
#> 
#> 3 variables in sample-level data:
#>   patient_id, neoadjuvant, overall_survival_days

As observed, by calling the tnbc object, information on the number of spots and genes per sample is displayed. For count statistics, the summarize_STlist function can be used:

summarize_STlist(tnbc)
#> # A tibble: 8 × 9
#>   sample_name spotscells genes min_cou…¹ mean_…² max_c…³ min_g…⁴ mean_…⁵ max_g…⁶
#>   <chr>            <int> <int>     <dbl>   <dbl>   <dbl>   <int>   <dbl>   <int>
#> 1 sample_117d       1949 23755       828  35264.  107737     533   6352.    9776
#> 2 sample_117e        754 21662      1022  34607.  121682     632   6141.   10091
#> 3 sample_120d       1359 21886       171  13831.   54618     150   3725.    7311
#> 4 sample_120e        844 20389        42   9997.   91663      41   3166.    9288
#> 5 sample_092a       1289 24227       559  43225.  103713     413   6951.    9650
#> 6 sample_092b       1376 23511      1413  34809.  117138     721   5959.    9680
#> 7 sample_094c       1542 24581        71  45805.  206674      70   7485.   11435
#> 8 sample_094d       1391 24167       298  37157.  130422     242   7061.   11129
#> # … with abbreviated variable names ¹​min_counts_per_spotcell,
#> #   ²​mean_counts_per_spotcell, ³​max_counts_per_spotcell,
#> #   ⁴​min_genes_per_spotcell, ⁵​mean_genes_per_spotcell, ⁶​max_genes_per_spotcell

The minimum number of counts per spot is 42, which seems low. We can look at the distribution of counts and genes per spot using the distribution_plots function:

cp <- distribution_plots(tnbc, plot_type='violin', plot_meta='total_counts')
cp[['total_counts']]

Now, let us remove spots with low counts by keeping only those spots with at least 5000 counts. We will also restrict the data set to spots expressing at least 1000 genes. We also will remove a few spots that have abnormally large number of counts, as in sample_094c. These criteria are not a rule, and samples in each study have to be carefully examined. For example, this criteria may be not enough to reduce the differences in counts, especially for sample_120d and samples_120e.

We run this filter with the filter_data function:

tnbc <- filter_data(tnbc, spot_minreads=5000, spot_mingenes=1000, spot_maxreads=150000)

cp2 <- distribution_plots(tnbc, plot_type='violin', plot_meta='total_counts')
cp2[['total_counts']]

Exploring variation between spatial arrays

The functions pseudobulk_samples and pseudobulk_pca_plot are the initial steps to obtain a quick snapshot of the variation in gene expression among samples. The function pseudobulk_samples creates (pseudo) “bulk” RNAseq data sets by combining all counts from each sample. Then, it log transforms the “pseudo bulk” RNAseq counts and runs a Principal Component Analysis (PCA). Note that the spatial coordinate information is not considered here, which is intended only as an exploratory analysis analysis. The max_var_genes argument is used to specify the maximum number of genes used for computation of PCA. The genes are selected based on their standard deviation across samples.

tnbc <- pseudobulk_samples(tnbc, max_var_genes=5000)

In this case, we apply the function to look for agreement between samples from the same patient: It is expected that tissue slices from the same patient are more similar among them than tissue slices from other patients. The pseudobulk_pca allows to map a sample-level variable to the points in the PCA by including the name of the column from the sample metadata (patient_id in this example).

pseudobulk_pca_plot(tnbc, plot_meta='patient_id')

Users can also generate a heatmap from pseudobulk counts by calling the pseudobulk_heatmap function, which also requires prior use of the pseudobulk_samples function. The number of variable genes to show can be controlled via the hm_display_genes argument.

hm_p <- pseudobulk_heatmap(tnbc, plot_meta='patient_id', hm_display_genes=30)

Transformation of spatially-resolved transcriptomics data

Many transformation methods are available for RNAseq count data. In spatialGE, the function transform_data applies log-transformation to the counts, after library size normalization performed on each sample separately. Similar to Seurat, it applies a scaling factor (scale_f=10000 by default).

tnbc <- transform_data(tnbc, method='log')

Users also can apply variance-stabilizing transformation (SCT; Hafemeister and Satija (2019)), which is another method increasingly used in single-cell and spatial transcriptomics studies. See here for details.

Visualization of gene expression from spatially-resolved transcriptomics data

After data transformation, expression of specific genes can be visualized using “quilt” plots. The function STplot shows the transformed expression of one or several genes. We have adopted the color palettes from the packages khroma and RColorBrewer. The name of a color palette can be passed using the argument color_pal. The default behavior of the function produces plots for all samples within the STlist, but we can pass specific samples to be plotted using the argument samples.

Let’s produce a quilt plot for the genes NDRG1 and IGKC (hypoxia and B-cell markers, respectively), for sample number 1 of patient 14 (samples=sample_094c).

quilts1 <- STplot(tnbc, 
                  genes=c('NDRG1', 'IGKC'), 
                  samples='sample_094c', 
                  color_pal='BuRd', 
                  ptsize=0.8)

Because spatialGE functions output lists of ggplot objects, we can plot the results side-by-side using functions such as ggarrange():

ggpubr::ggarrange(plotlist=quilts1, nrow=1, ncol=2, legend='bottom')

We can see that gene expression patterns of both genes are non-overlapping: IGKC is expressed in the upper right portion of the tissue, whereas NDRG1 is expressed in the right bottom portion (although relatively lower compared to IGKC as indicated by the white-colored spots). The location of gene expression of those two genes may be indicative of an immune-infiltrated area and a tumor area. With the help of spatial interpolation, visualization of these regions can be easier as will be showed next.

Spatial interpolation of gene expression

We can predict a smooth gene expression surface for each sample. In spatialGE, this prediction is achieved by using a spatial interpolation method very popular in spatial statistics. The method known as ‘kriging’ allows the estimation of gene expression values in the un-sampled areas between spots, or cells/spots that were filtered during data quality control. Estimating a transcriptomic surface via kriging assumes that gene expression of two given points is correlated to the spatial distance between them.

The function gene_interpolation performs kriging of gene expression via the fields package. We specify that kriging will be performed for two of the spatial samples (samples=c('sample_094c', 'sample_117e')):

tnbc <- gene_interpolation(tnbc, 
                           genes=c('NDRG1', 'IGKC'),
                           samples=c('sample_094c', 'sample_117e'))
#> Gene interpolation started.
#> Gene interpolation completed in 1.59 min.
#> 

Generating gene expression surfaces can be time consuming. The finer resolution to which the surface is to be predicted (ngrid argument), the longer the time it takes. The execution time also depends on the number of spots/cells. The surfaces can be visualized using the STplot_interpolation() function:

kriges1 <- STplot_interpolation(tnbc,
                                genes=c('NDRG1', 'IGKC'),
                                samples='sample_094c')
#> Loading required namespace: rgeos
#> Loading required package: ggpolypath
#> Loading required package: ggplot2
ggpubr::ggarrange(plotlist=kriges1, nrow=1, ncol=2, common.legend=T, legend='bottom')

By looking at the transcriptomic surfaces of the two genes, it is easier to detect where “pockets” of high and low expression are located within the tissue. It is now more evident that expression of NDRG1 is higher in the lower right region of the tumor slice, as well in a smaller area to the left (potentially another hypoxic region).

If tissue images were uploaded to the STlist (tissue_hires_image.png.gz), it can be displayed for comparison using the plot_image function:

ip = plot_image(tnbc, samples='sample_094c')

ip[['image_sample_094c']]

Unsupervised spatially-informed clustering (STclust)

Detecting tissue compartments or niches is an important part of the study of the tissue architecture. We can do this by applying STclust, a spatially-informed clustering method implemented in spatialGE. The STclust method uses weighted average matrices to capture the transcriptomic differences among the cells/spots. As a first step in STclust, top variable genes are identified via Seurat’s FindVariableFeatures, and transcriptomic scaled distances are calculated using only those genes. Next, scaled euclidean distances are computed from the spatial coordinates of the spots/cells. The user defines a weight (ws) from 0 to 1, to apply to the physical distances. The higher the weight, the less biologically meaningful the clustering solution is, given that the clusters would only reflect the physical distances between the spots/cells and less information on the transcriptomic profiles will be used. After many tests, we have found that weights between 0.025 - 0.25 are enough to capture the tissue heterogeneity. By default, STclust uses dynamic tree cuts (Langfelder, Zhang, and Horvath 2007) to define the number of clusters. But users can also test a series of k values (ks). For a more detailed description of the method, please refer to the paper describing spatialGE and STclust (Ospina et al. 2022).

We’ll try several weights to see it’s effect on the cluster assignments:

tnbc <- STclust(tnbc, 
                ks='dtc', 
                ws=c(0, 0.025, 0.05, 0.2))
#> STclust started.
#> STclust completed in 0.85 min.
#> 

Results of clustering can be plotted via the STplot function:

cluster_p <- STplot(tnbc, 
                    samples='sample_094c', 
                    ws=c(0, 0.025, 0.05, 0.2),
                    color_pal='highcontrast')
ggpubr::ggarrange(plotlist=cluster_p, nrow=2, ncol=2, legend='right')

We can see that from w=0 and w=0.05, we can only detect two tissue niches. At w=0.025, we gain higher resolution as one of the clusters is split and a third (‘yellow’) cluster appears, potentially indicating that an different transcriptional profile is present there. At w=0.2, the clusters seem too compact, indicating that the weight of spatial information is probably too high. We have used here the dynamic tree cuts (dtc) to automatically select the number of clusters, resulting in very coarse resolution tissue niches, however, users can define their own range of k to be evaluated, allowing further detection of tissue compartments.

Association between spatial heterogeneity and sample-level variables

To explore the relationship between a clinical (sample-level) variable of interest and the level of gene expression spatial uniformity within a sample, we can use the SThet() function:

tnbc <- SThet(tnbc, 
              genes=c('NDRG1', 'IGKC'),
              method='moran')
#> SThet started.
#>    Calculating spatial weights...
#> SThet completed in 0.21 min.
#> 

The SThet function calculates the Moran’s I statsitic (or Geary’s C) to measure the level of spatial heterogeneity in the expression of the genes ( NDRG1, IGKC). The estimates can be compared across samples using the function compare_SThet()

p <- compare_SThet(tnbc, 
                   samplemeta='overall_survival_days', 
                   color_by='patient_id',
                   gene=c('NDRG1', 'IGKC'), 
                   color_pal='muted',
                   ptsize=3)

p

The calculation of spatial statistics with SThet and and multi-sample comparison with compare_SThet provides and easy way to identify samples and genes exhibiting spatial patterns. The previous figure shows that expression of NDRG1 is more spatially uniform (lower Moran’s I) across the tissues in samples from patients 2 and 8 compared to patents 9 and 14. The samples with higher spatial uniformity in the expression of NDRG1 also tended to have higher overall survival. Trends are less clear for IGKC, however, it looks like samples where expression of IGKC was spatially aggregated in “pockets” (higher Moran’s I) tended to have lower survival (but note that patient 9 does not follow this trend). As studies using spatial transcriptomics become larger, more samples will provide more insightful patterns into the association of gene expression spatial distribution and non-spatial traits associated with the tissues.

The computed statistics are stored in the STlist for additional analysis/plotting that the user may want to complete. The statistics value can be accessed as a data frame using the get_gene_meta function:

get_gene_meta(tnbc, sthet_only=T)
#> # A tibble: 16 × 10
#>    sample  patie…¹ neoad…² overa…³ gene  gene_…⁴ gene_…⁵ vst.v…⁶ moran_i geary_c
#>    <chr>     <dbl> <chr>     <dbl> <chr>   <dbl>   <dbl>   <dbl>   <dbl> <lgl>  
#>  1 sample…       2 yes        2733 IGKC    3.07    0.787    9.73  0.0778 NA     
#>  2 sample…       2 yes        2733 NDRG1   2.13    0.436    1.74  0.0818 NA     
#>  3 sample…       2 yes        2733 IGKC    2.54    0.799   12.6   0.0813 NA     
#>  4 sample…       2 yes        2733 NDRG1   2.45    0.502    1.66  0.157  NA     
#>  5 sample…       8 no         1453 IGKC    2.98    1.37    14.0   0.133  NA     
#>  6 sample…       8 no         1453 NDRG1   0.672   0.514    3.11  0.0303 NA     
#>  7 sample…       8 no         1453 IGKC    2.77    1.23     7.80  0.119  NA     
#>  8 sample…       8 no         1453 NDRG1   1.05    0.714    1.97  0.131  NA     
#>  9 sample…       9 no          308 IGKC    0.371   0.626   11.4   0.0916 NA     
#> 10 sample…       9 no          308 NDRG1   1.11    0.734    5.89  0.237  NA     
#> 11 sample…       9 no          308 IGKC    0.105   0.412    9.54  0.0561 NA     
#> 12 sample…       9 no          308 NDRG1   1.63    0.872    4.08  0.116  NA     
#> 13 sample…      14 yes         718 IGKC    4.06    1.60     6.31  0.203  NA     
#> 14 sample…      14 yes         718 NDRG1   1.38    0.780    6.04  0.199  NA     
#> 15 sample…      14 yes         718 IGKC    3.23    1.43     9.12  0.164  NA     
#> 16 sample…      14 yes         718 NDRG1   1.25    0.768    6.69  0.174  NA     
#> # … with abbreviated variable names ¹​patient_id, ²​neoadjuvant,
#> #   ³​overall_survival_days, ⁴​gene_mean, ⁵​gene_stdevs,
#> #   ⁶​vst.variance.standardized
How can the statistics generated by SThet can be interpreted?

See the table below for a simplistic interpretation of the spatial autocorrelation statistics calculated in spatialGE:

Statistic Clustered expression No expression pattern Uniform expression
Moran’s I Closer to 1 Closer to 0 Closer to -1
Geary’s C Closer to 0 Closer to 1 Closer to 2
Note: The boundaries indicated for each statistic are reached when the number of spots/cells is very large.

To better understand how the Moran’s I and Geary’s C statistics quantify spatial heterogeneity, tissue can be simulated using the sc.SpatialSIM and spatstat packages. Also tidyverse and janitor for some data manipulation.

The sc.SpatialSIM package uses spatial point processes to simulate the locations of spots/cells within a tissue. To facilitate interpretability of the Moran’s I and Geary’s C statistics, first tissue will be simulated, and then gene expression values will be simulated.

The first step is to create a spatstat observation window:

wdw <- owin(xrange=c(0, 3), yrange=c(0, 3))
sim_visium <- CreateSimulationObject(sims=1, cell_types=1, window=wdw)

Next, the sc.SpatialSIM is used to generate the spatial point process (positions of the spots) within the observation window. Then, assignments of spots to tissue domains are simulated and visualized:

# Generate point process
# Then, simulate tissue compartments
set.seed(12345)
sim_visium <- GenerateSpatialPattern(sim_visium, gridded=T, overwrite=T) %>% 
  GenerateTissue(k=1, xmin=1, xmax=2, ymin=1, ymax=2, sdmin=1, sdmax=2, overwrite=T)
#> x and y range inside window boundary
#> Computing tissue probability

PlotSimulation(sim_visium, which=1, what="tissue points") + 
  scale_shape_manual(values=c(19, 1))

Now, the simulated tissue domain assignments are extracted from the SpatSimObj object. Gene counts will be simulated in such a way that:

  • Expression of “gene_1” is exclusive to “Tissue_1”. If a spot was assigned to “Tissue 1”, expression of “gene_1” is 1 (high expression). If assigned to “Tissue 2”, expression of “gene_1” is 0.1 (low expression). Given “Tissue 1” spots are aggregated towards the center of the tissue, it is expected that Moran’s I > 0 and Geary’s C < 1 for “gene_1”.
  • Spots with high expression of “gene_2” (expression = 1) are equally separated from spots with low expression (expression = 0.1). This pattern is highly unlikely in a biological tissue. It is expected this pattern yields a Moran’s I < 0 and a Geary’s > 1.
  • Expression of “gene_3” results from randomly assigning expression of “gene_2” across the entire tissue. This pattern should result in Moran’s I closer to 0 and Geary’s C closer to 1.
# Extract tissue assignments from the `SpatSimObj` object
# Simulate expression of 'gene_1'
sim_visium_df <- sim_visium@`Spatial Files`[[1]] %>%
  clean_names() %>% 
  mutate(gene_1=case_when(tissue_assignment == 'Tissue 1' ~ 1, TRUE ~ 0.1))

# Generate expression patter of "gene_2"
for(i in 1:nrow(sim_visium_df)){
  if(i%%2 == 0){
    sim_visium_df[i, 'gene_2'] = 1
  } else{
    sim_visium_df[i, 'gene_2'] = 0.1
  }
}

# Generate expression of "gene_3"
# Set seed for resproducibility
set.seed(12345)
sim_visium_df[['gene_3']] <- sample(sim_visium_df[['gene_2']])

To visualize the simulated expression and run SThet, an STlist is created:

# Extract simulated expression data
sim_expr <- sim_visium_df %>% 
  add_column(libname=paste0('spot', seq(1:nrow(.)))) %>%
  select(c('libname', 'gene_1', 'gene_2', 'gene_3')) %>%
  column_to_rownames('libname') %>% t() %>% 
  as.data.frame() %>% rownames_to_column('genename')

# Extract simulated spot locations
sim_xy <- sim_visium_df %>% 
  add_column(libname=paste0('spot', seq(1:nrow(.)))) %>%
  select(c('libname', 'y', 'x'))

# The `STlist` function can take a list of data frames
simulated <- STlist(rnacounts=list(sim_visium=sim_expr), 
                   spotcoords=list(sim_visium=sim_xy))
#> Found list of dataframes.
#> Matching gene expression and coordinate data...
#> Converting counts to sparse matrices
#> Completed STlist!
#> 

# Plot expression
ps <- STplot(simulated, genes=c('gene_1', 'gene_2', 'gene_3'), data_type='raw', color_pal='sunset', ptsize=1)
ggpubr::ggarrange(plotlist=ps, ncol=3)

As expected, simulated expression of “gene_1” is aggregated (i.e., clustered, red spots concentrated towards the center of the tissue). This result was intently obtained by setting the simulation parameters to xmin=2, xmax=3, ymin=2, ymax=3 in sc.SpatialSIM.

# The `SThet` function requires normalized data
simulated <- transform_data(simulated)

# Run `SThet`
simulated <- SThet(simulated, genes=c('gene_1', 'gene_2', 'gene_3'), method=c('moran', 'geary'))
#> SThet started.
#>    Calculating spatial weights...
#> SThet completed in 0.02 min.
#> 

# Extract results
get_gene_meta(simulated, sthet_only=T)
#> # A tibble: 3 × 6
#>   sample     gene   gene_mean gene_stdevs  moran_i geary_c
#>   <chr>      <chr>      <dbl>       <dbl>    <dbl>   <dbl>
#> 1 sim_visium gene_1      7.80       0.974  0.0193    0.966
#> 2 sim_visium gene_2      7.68       0.957 -0.0156    1.01 
#> 3 sim_visium gene_3      7.68       1.01   0.00371   0.987

The results for each of the metrics are as expected: “gene_1” shows aggregation/clustering, indicative of “hot-spot” expression. “gene_2” and “gene_3” show uniform and random expression respectively. Notice that none of these values are very far from the “random” expectation (Moran’s I = 0 and Geary’s C = 1). One of reasons for this result is the effect of library-size normalization. In addition, obtaining extreme values of I and C require extreme spatial patterns, which are unlikely to be observed even in the way data has been simulated here.

References

Bassiouni, Rania, Michael O. Idowu, Lee D. Gibbs, Valentina Robila, Pamela J. Grizzard, Michelle G. Webb, Jiarong Song, Ashley Noriega, David W. Craig, and John D. Carpten. 2023. “Spatial Transcriptomic Analysis of a Diverse Patient Cohort Reveals a Conserved Architecture in Triple-Negative Breast Cancer.” Cancer Research 83: 34–48. https://doi.org/10.1158/0008-5472.CAN-22-2682.
Hafemeister, Christoph, and Rahul Satija. 2019. “Normalization and Variance Stabilization of Single-Cell RNA-Seq Data Using Regularized Negative Binomial Regression.” Genome Biology 20: 296. https://doi.org/10.1186/s13059-019-1874-1.
Langfelder, Peter, Bin Zhang, and Steve Horvath. 2007. Defining clusters from a hierarchical cluster tree: the Dynamic Tree Cut package for R.” Bioinformatics 24: 719–20. https://doi.org/10.1093/bioinformatics/btm563.
Ospina, Oscar E, Christopher M Wilson, Alex C Soupir, Anders Berglund, Inna Smalley, Kenneth Y Tsai, and Brooke L Fridley. 2022. spatialGE: quantification and visualization of the tumor microenvironment heterogeneity using spatial transcriptomics.” Bioinformatics 38: 2645–47. https://doi.org/10.1093/bioinformatics/btac145.
Session Info
sessionInfo()
#> R version 4.1.2 (2021-11-01)
#> Platform: x86_64-apple-darwin17.0 (64-bit)
#> Running under: macOS Big Sur 10.16
#> 
#> Matrix products: default
#> BLAS:   /Library/Frameworks/R.framework/Versions/4.1/Resources/lib/libRblas.0.dylib
#> LAPACK: /Library/Frameworks/R.framework/Versions/4.1/Resources/lib/libRlapack.dylib
#> 
#> locale:
#> [1] en_US.UTF-8/en_US.UTF-8/en_US.UTF-8/C/en_US.UTF-8/en_US.UTF-8
#> 
#> attached base packages:
#> [1] stats     graphics  grDevices utils     datasets  methods   base     
#> 
#> other attached packages:
#>  [1] janitor_2.1.0            forcats_0.5.1            stringr_1.5.0           
#>  [4] dplyr_1.1.0              purrr_1.0.1              readr_2.1.2             
#>  [7] tidyr_1.3.0              tibble_3.2.1             tidyverse_1.3.1         
#> [10] sc.SpatialSIM_0.0.0.9000 spatstat_2.3-0           spatstat.linnet_2.3-1   
#> [13] spatstat.core_2.3-2      rpart_4.1-15             nlme_3.1-153            
#> [16] spatstat.geom_3.0-6      spatstat.data_3.0-0      magrittr_2.0.2          
#> [19] ggpolypath_0.1.0         ggplot2_3.3.5            fields_14.1             
#> [22] viridis_0.6.2            viridisLite_0.4.0        spam_2.9-1              
#> [25] spatialGE_1.2.0         
#> 
#> loaded via a namespace (and not attached):
#>   [1] scattermore_0.7       SeuratObject_4.0.4    ragg_1.2.5           
#>   [4] bit64_4.0.5           knitr_1.37            irlba_2.3.5          
#>   [7] data.table_1.14.2     doParallel_1.0.16     generics_0.1.2       
#>  [10] BiocGenerics_0.40.0   terra_1.7-3           cowplot_1.1.1        
#>  [13] RANN_2.6.1            proxy_0.4-26          future_1.23.0        
#>  [16] wordspace_0.2-6       bit_4.0.4             tzdb_0.3.0           
#>  [19] wk_0.6.0              webshot_0.5.2         xml2_1.3.3           
#>  [22] lubridate_1.8.0       httpuv_1.6.5          assertthat_0.2.1     
#>  [25] xfun_0.31             hms_1.1.1             jquerylib_0.1.4      
#>  [28] evaluate_0.14         promises_1.2.0.1      fansi_1.0.2          
#>  [31] dbplyr_2.1.1          readxl_1.3.1          igraph_1.2.11        
#>  [34] DBI_1.1.2             htmlwidgets_1.5.4     sparsesvd_0.2        
#>  [37] reshape_0.8.9         stats4_4.1.2          spdep_1.2-2          
#>  [40] ellipsis_0.3.2        ggpubr_0.4.0          backports_1.4.1      
#>  [43] V8_4.0.0              deldir_1.0-6          vctrs_0.6.2          
#>  [46] Cairo_1.6-0           ROCR_1.0-11           abind_1.4-5          
#>  [49] cachem_1.0.6          withr_2.5.0           ggforce_0.3.3        
#>  [52] rgdal_1.6-4           vroom_1.5.7           sctransform_0.3.3    
#>  [55] goftest_1.2-3         svglite_2.0.0         cluster_2.1.2        
#>  [58] dotCall64_1.0-1       lazyeval_0.2.2        crayon_1.4.2         
#>  [61] hdf5r_1.3.8           units_0.7-2           pkgconfig_2.0.3      
#>  [64] labeling_0.4.2        tweenr_1.0.2          rlang_1.1.1          
#>  [67] globals_0.14.0        lifecycle_1.0.3       miniUI_0.1.1.1       
#>  [70] modelr_0.1.8          cellranger_1.1.0      rprojroot_2.0.2      
#>  [73] polyclip_1.10-0       matrixStats_0.61.0    lmtest_0.9-39        
#>  [76] Matrix_1.4-1          raster_3.5-15         carData_3.0-5        
#>  [79] boot_1.3-28           zoo_1.8-9             reprex_2.0.1         
#>  [82] ggridges_0.5.3        GlobalOptions_0.1.2   png_0.1-7            
#>  [85] rjson_0.2.21          KernSmooth_2.23-20    shape_1.4.6          
#>  [88] classInt_0.4-3        s2_1.1.2              parallelly_1.30.0    
#>  [91] iotools_0.3-2         rstatix_0.7.0         S4Vectors_0.32.3     
#>  [94] ggsignif_0.6.3        scales_1.2.0          memoise_2.0.1        
#>  [97] plyr_1.8.6            ica_1.0-2             compiler_4.1.2       
#> [100] kableExtra_1.3.4      concaveman_1.1.0      RColorBrewer_1.1-2   
#> [103] clue_0.3-60           fitdistrplus_1.1-6    snakecase_0.11.0     
#> [106] cli_3.6.0             listenv_0.8.0         patchwork_1.1.1      
#> [109] pbapply_1.5-0         MASS_7.3-54           mgcv_1.8-38          
#> [112] tidyselect_1.2.0      stringi_1.7.6         textshaping_0.3.6    
#> [115] highr_0.9             yaml_2.2.2            ggrepel_0.9.1        
#> [118] pbmcapply_1.5.1       grid_4.1.2            sass_0.4.0           
#> [121] tools_4.1.2           future.apply_1.8.1    parallel_4.1.2       
#> [124] circlize_0.4.14       rstudioapi_0.13       khroma_1.8.0         
#> [127] foreach_1.5.2         gridExtra_2.3         farver_2.1.0         
#> [130] Rtsne_0.15            digest_0.6.29         rgeos_0.5-9          
#> [133] shiny_1.7.1           Rcpp_1.0.10           car_3.0-12           
#> [136] broom_0.7.12          later_1.3.0           RcppAnnoy_0.0.19     
#> [139] httr_1.4.2            sf_1.0-9              ComplexHeatmap_2.11.1
#> [142] colorspace_2.0-2      rvest_1.0.2           fs_1.5.2             
#> [145] tensor_1.5            reticulate_1.24       IRanges_2.28.0       
#> [148] splines_4.1.2         uwot_0.1.11           spatstat.utils_3.0-1 
#> [151] pkgdown_2.0.7         sp_1.4-6              plotly_4.10.0        
#> [154] spData_2.0.1          systemfonts_1.0.4     xtable_1.8-4         
#> [157] jsonlite_1.7.3        dynamicTreeCut_1.63-1 R6_2.5.1             
#> [160] pillar_1.8.1          htmltools_0.5.2       mime_0.12            
#> [163] glue_1.6.2            fastmap_1.1.0         class_7.3-19         
#> [166] codetools_0.2-18      maps_3.4.0            utf8_1.2.2           
#> [169] lattice_0.20-45       bslib_0.3.1           spatstat.sparse_3.0-0
#> [172] curl_4.3.2            leiden_0.3.9          magick_2.7.3         
#> [175] survival_3.4-0        rmarkdown_2.14        desc_1.4.2           
#> [178] munsell_0.5.0         e1071_1.7-9           GetoptLong_1.0.5     
#> [181] iterators_1.0.13      haven_2.4.3           reshape2_1.4.4       
#> [184] gtable_0.3.0          Seurat_4.1.0