Title: | R wrapper for SpatialDE |
---|---|
Description: | SpatialDE is a method to find spatially variable genes (SVG) from spatial transcriptomics data. This package provides wrappers to use the Python SpatialDE library in R, using reticulate and basilisk. |
Authors: | Davide Corso [aut] , Milan Malfait [aut] , Lambda Moses [aut] , Gabriele Sales [cre] |
Maintainer: | Gabriele Sales <[email protected]> |
License: | MIT + file LICENSE |
Version: | 1.13.0 |
Built: | 2024-10-31 05:29:00 UTC |
Source: | https://github.com/bioc/spatialDE |
Plot Fraction Spatial Variance vs Q-value
FSV_sig( results, ms_results = NULL, certain_only = FALSE, log_x = FALSE, do_label = TRUE, covariate_names = NULL )
FSV_sig( results, ms_results = NULL, certain_only = FALSE, log_x = FALSE, do_label = TRUE, covariate_names = NULL )
results |
results from SpatialDE. |
ms_results |
model selection results, should be a data frame with
columns |
certain_only |
only plot results with narrow 95% confidence interval. |
log_x |
Whether to display x axis in log scale. |
do_label |
display gene names for statistically significant genes,
default |
covariate_names |
names of covariates as a reference, default to |
A ggplot2
object.
Davide Corso, Milan Malfait, Lambda Moses
Svensson, V., Teichmann, S. & Stegle, O. SpatialDE: identification of spatially variable genes. Nat Methods 15, 343–346 (2018). https://doi.org/10.1038/nmeth.4636
SpatialDE 1.1.3: the version of the Python package used under the hood.
## Mock up a SpatialExperiment object wit 400 cells and 3 genes set.seed(42) spe <- mockSVG(size = 20, tot_genes = 3, de_genes = 1, return_SPE = TRUE) ## Run spatialDE with S4 integration results <- spatialDE(spe) ## Run model search msearch <- modelSearch(spe, de_results = results, qval_thresh = NULL, verbose = FALSE) plot <- FSV_sig(results, msearch)
## Mock up a SpatialExperiment object wit 400 cells and 3 genes set.seed(42) spe <- mockSVG(size = 20, tot_genes = 3, de_genes = 1, return_SPE = TRUE) ## Run spatialDE with S4 integration results <- spatialDE(spe) ## Run model search msearch <- modelSearch(spe, de_results = results, qval_thresh = NULL, verbose = FALSE) plot <- FSV_sig(results, msearch)
Coordinates and total counts for the samples from the Mouse Olfactory Bulb data generated by Stahl et al. (2016). This data was originally downloaded from https://github.com/Teichlab/SpatialDE/blob/master/Analysis/MouseOB/MOB_sample_info.csv.
data(MOB_sample_info)
data(MOB_sample_info)
A data.frame
with 262 rows and 3 variables as columns: the x
and
y
coordinates and total_counts
corresponding to each spot.
Ståhl, P. L. et al. (2016) 'Visualization and analysis of gene expression in tissue sections by spatial transcriptomics', Science, 353(6294), p. 78. doi: 10.1126/science.aaf2403.
Generate count matrix for spatially variable genes.
mockSVG(size, tot_genes, de_genes, return_SPE = FALSE)
mockSVG(size, tot_genes, de_genes, return_SPE = FALSE)
size |
An |
tot_genes |
An |
de_genes |
An |
return_SPE |
A |
If return_SPE = TRUE
, returns a SpatialExperiment object.
If not, a list
containing:
coordinates
: data.frame
with x
and y
columns;
counts
: matrix
with generated gene counts.
spe <- mockSVG(size = 20, tot_genes = 3, de_genes = 1, return_SPE = TRUE) spe
spe <- mockSVG(size = 20, tot_genes = 3, de_genes = 1, return_SPE = TRUE) spe
Classify DE genes to interpretable fitting classes.
model_search(x, coordinates, de_results, qval_thresh = 0.05, verbose = FALSE)
model_search(x, coordinates, de_results, qval_thresh = 0.05, verbose = FALSE)
x |
|
coordinates |
|
de_results |
|
qval_thresh |
|
verbose |
|
data.frame
of model_search results.
Svensson, V., Teichmann, S. & Stegle, O. SpatialDE: identification of spatially variable genes. Nat Methods 15, 343–346 (2018). https://doi.org/10.1038/nmeth.4636
## Mock up a SpatialExperiment object wit 400 cells and 3 genes set.seed(42) mock <- mockSVG(size = 20, tot_genes = 3, de_genes = 1) stabilized <- stabilize(mock$counts) sample_info <- mock$coordinates sample_info$total_counts <- colSums(mock$counts) regressed <- regress_out(counts = stabilized, sample_info = sample_info) ## Run SpatialDE de_results <- run(regressed, coordinates = mock$coordinates) ## Run model search ms_results <- model_search( x = regressed, coordinates = mock$coordinates, de_results = de_results, qval_thresh = NULL )
## Mock up a SpatialExperiment object wit 400 cells and 3 genes set.seed(42) mock <- mockSVG(size = 20, tot_genes = 3, de_genes = 1) stabilized <- stabilize(mock$counts) sample_info <- mock$coordinates sample_info$total_counts <- colSums(mock$counts) regressed <- regress_out(counts = stabilized, sample_info = sample_info) ## Run SpatialDE de_results <- run(regressed, coordinates = mock$coordinates) ## Run model search ms_results <- model_search( x = regressed, coordinates = mock$coordinates, de_results = de_results, qval_thresh = NULL )
Compare model fits with different models, using the SpatialDE Python package.
modelSearch(x, de_results, ...) ## S4 method for signature 'matrix' modelSearch(x, de_results, coordinates, qval_thresh = 0.05, verbose = FALSE) ## S4 method for signature 'SpatialExperiment' modelSearch( x, de_results, assay_type = "counts", qval_thresh = 0.05, verbose = FALSE )
modelSearch(x, de_results, ...) ## S4 method for signature 'matrix' modelSearch(x, de_results, coordinates, qval_thresh = 0.05, verbose = FALSE) ## S4 method for signature 'SpatialExperiment' modelSearch( x, de_results, assay_type = "counts", qval_thresh = 0.05, verbose = FALSE )
x |
A numeric Alternatively, a SpatialExperiment object. |
de_results |
|
... |
For the generic, arguments to pass to specific methods. |
coordinates |
A For the SpatialExperiment method, coordinates are taken from
|
qval_thresh |
|
verbose |
A |
assay_type |
A |
data.frame
of model_search results.
Davide Corso, Milan Malfait, Lambda Moses
Svensson, V., Teichmann, S. & Stegle, O. SpatialDE: identification of spatially variable genes. Nat Methods 15, 343–346 (2018). https://doi.org/10.1038/nmeth.4636
SpatialDE 1.1.3: the version of the Python package used under the hood.
The individual steps performed by this function: stabilize()
,
regress_out()
and model_search()
.
## Mock up a SpatialExperiment object wit 400 cells and 3 genes set.seed(42) spe <- mockSVG(size = 20, tot_genes = 3, de_genes = 1, return_SPE = TRUE) ## Run spatialDE with S4 integration de_results <- spatialDE(spe) ## Run model search model_search <- modelSearch(spe, de_results = de_results, qval_thresh = NULL, verbose = FALSE )
## Mock up a SpatialExperiment object wit 400 cells and 3 genes set.seed(42) spe <- mockSVG(size = 20, tot_genes = 3, de_genes = 1, return_SPE = TRUE) ## Run spatialDE with S4 integration de_results <- spatialDE(spe) ## Run model search model_search <- modelSearch(spe, de_results = de_results, qval_thresh = NULL, verbose = FALSE )
Plot Spatial Patterns of Multiple Genes
multiGenePlots(x, ...) ## S4 method for signature 'matrix' multiGenePlots( x, coordinates, genes_plot, viridis_option = "D", ncol = 2, point_size = 1, dark_theme = TRUE ) ## S4 method for signature 'SpatialExperiment' multiGenePlots( x, assay_type = "counts", genes_plot, viridis_option = "D", ncol = 2, point_size = 1, dark_theme = TRUE )
multiGenePlots(x, ...) ## S4 method for signature 'matrix' multiGenePlots( x, coordinates, genes_plot, viridis_option = "D", ncol = 2, point_size = 1, dark_theme = TRUE ) ## S4 method for signature 'SpatialExperiment' multiGenePlots( x, assay_type = "counts", genes_plot, viridis_option = "D", ncol = 2, point_size = 1, dark_theme = TRUE )
x |
A numeric Alternatively, a SpatialExperiment object. |
... |
For the generic, arguments to pass to specific methods. |
coordinates |
A For the SpatialExperiment method, coordinates are taken from
|
genes_plot |
character vector specifying which genes are to be plotted. |
viridis_option |
This function uses the |
ncol |
Number of columns to arrange the plots. |
point_size |
Point size of each plot. |
dark_theme |
Whether dark background should be used; this is helpful to
highlight cells with high expression when using the |
assay_type |
A |
This function draws a plot for each specified genes
Davide Corso, Milan Malfait, Lambda Moses
Svensson, V., Teichmann, S. & Stegle, O. SpatialDE: identification of spatially variable genes. Nat Methods 15, 343–346 (2018). https://doi.org/10.1038/nmeth.4636
SpatialDE 1.1.3: the version of the Python package used under the hood.
The individual steps performed by this function: stabilize()
,
spatialDE()
.
For further analysis of the DE results:
model_search()
and spatial_patterns()
.
## Mock up a SpatialExperiment object wit 400 cells and 3 genes set.seed(42) spe <- mockSVG(size = 20, tot_genes = 3, de_genes = 1, return_SPE = TRUE) ## Run spatialDE results <- spatialDE(spe) ordered_spe_results <- results[order(results$qval), ] head(ordered_spe_results) plots <- multiGenePlots(spe, assay_type = "counts", ordered_spe_results$g, point_size = 4, viridis_option = "D" )
## Mock up a SpatialExperiment object wit 400 cells and 3 genes set.seed(42) spe <- mockSVG(size = 20, tot_genes = 3, de_genes = 1, return_SPE = TRUE) ## Run spatialDE results <- spatialDE(spe) ordered_spe_results <- results[order(results$qval), ] head(ordered_spe_results) plots <- multiGenePlots(spe, assay_type = "counts", ordered_spe_results$g, point_size = 4, viridis_option = "D" )
Regresses out the effect of library size.
This function is a wrapper for regress_out
from the
NaiveDE Python package.
regress_out(counts, sample_info)
regress_out(counts, sample_info)
counts |
|
sample_info |
|
matrix
of normalized counts.
## Mock up a SpatialExperiment object wit 400 cells and 3 genes set.seed(42) mock <- mockSVG(20, 3, 1) stabilized <- stabilize(mock$counts) sample_info <- mock$coordinates sample_info$total_counts <- colSums(mock$counts) regressed <- regress_out(counts = stabilized, sample_info = sample_info)
## Mock up a SpatialExperiment object wit 400 cells and 3 genes set.seed(42) mock <- mockSVG(20, 3, 1) stabilized <- stabilize(mock$counts) sample_info <- mock$coordinates sample_info$total_counts <- colSums(mock$counts) regressed <- regress_out(counts = stabilized, sample_info = sample_info)
Replicate 11 from the spatially dependent gene expression data from the mouse olfactory bulb generated by Stahl et al. (2016). This data was originally downloaded from https://github.com/Teichlab/SpatialDE/blob/master/Analysis/MouseOB/data/Rep11_MOB_0.csv.
data(Rep11_MOB_0)
data(Rep11_MOB_0)
A matrix with 16218 genes as rows and 262 spots as columns.
Ståhl, P. L. et al. (2016) 'Visualization and analysis of gene expression in tissue sections by spatial transcriptomics', Science, 353(6294), p. 78. doi: 10.1126/science.aaf2403.
Wraps the run
function from the
SpatialDE Python package.
run(x, coordinates, verbose = FALSE)
run(x, coordinates, verbose = FALSE)
x |
|
coordinates |
|
verbose |
|
A data.frame
with DE results where each row is a gene and columns
contain relevant statistics.
The most important columns are:
g
: the name of the gene
pval
: the p-value for spatial differential expression
qval
: the q-value, indicating significance after correcting for
multiple testing
l
: A parameter indicating the distance scale a gene changes expression
over
Svensson, V., Teichmann, S. & Stegle, O. SpatialDE: identification of spatially variable genes. Nat Methods 15, 343–346 (2018). https://doi.org/10.1038/nmeth.4636
## Mock up a SpatialExperiment object wit 400 cells and 3 genes set.seed(42) mock <- mockSVG(size = 20, tot_genes = 3, de_genes = 1) stabilized <- stabilize(mock$counts) sample_info <- mock$coordinates sample_info$total_counts <- colSums(mock$counts) regressed <- regress_out(counts = stabilized, sample_info = sample_info) ## Run SpatialDE de_results <- run(regressed, coordinates = mock$coordinates)
## Mock up a SpatialExperiment object wit 400 cells and 3 genes set.seed(42) mock <- mockSVG(size = 20, tot_genes = 3, de_genes = 1) stabilized <- stabilize(mock$counts) sample_info <- mock$coordinates sample_info$total_counts <- colSums(mock$counts) regressed <- regress_out(counts = stabilized, sample_info = sample_info) ## Run SpatialDE de_results <- run(regressed, coordinates = mock$coordinates)
Group spatially variable genes into spatial patterns using automatic expression histology (AEH)
spatial_patterns( x, coordinates, de_results, qval_thresh = 0.05, n_patterns, length, verbose = FALSE )
spatial_patterns( x, coordinates, de_results, qval_thresh = 0.05, n_patterns, length, verbose = FALSE )
x |
|
coordinates |
|
de_results |
|
qval_thresh |
|
n_patterns |
|
length |
|
verbose |
|
list
of two dataframe (pattern_results, patterns):
pattern_results
dataframe with pattern membership information for each
gene.
patterns
the posterior mean underlying expression fro genes in given
spatial patterns.
Svensson, V., Teichmann, S. & Stegle, O. SpatialDE: identification of spatially variable genes. Nat Methods 15, 343–346 (2018). https://doi.org/10.1038/nmeth.4636
## Mock up a SpatialExperiment object wit 400 cells and 3 genes set.seed(42) mock <- mockSVG(size = 20, tot_genes = 3, de_genes = 1) stabilized <- stabilize(mock$counts) sample_info <- mock$coordinates sample_info$total_counts <- colSums(mock$counts) regressed <- regress_out(counts = stabilized, sample_info = sample_info) ## Run SpatialDE de_results <- run(x = regressed, coordinates = mock$coordinates) ## Run Spatial_patterns sp <- spatial_patterns( x = regressed, coordinates = mock$coordinates, de_results = de_results, qval_thresh = NULL, n_patterns = 5, length = 1.5 ) sp$pattern_results sp$patterns
## Mock up a SpatialExperiment object wit 400 cells and 3 genes set.seed(42) mock <- mockSVG(size = 20, tot_genes = 3, de_genes = 1) stabilized <- stabilize(mock$counts) sample_info <- mock$coordinates sample_info$total_counts <- colSums(mock$counts) regressed <- regress_out(counts = stabilized, sample_info = sample_info) ## Run SpatialDE de_results <- run(x = regressed, coordinates = mock$coordinates) ## Run Spatial_patterns sp <- spatial_patterns( x = regressed, coordinates = mock$coordinates, de_results = de_results, qval_thresh = NULL, n_patterns = 5, length = 1.5 ) sp$pattern_results sp$patterns
Identify genes that significantly depend on spatial coordinates with the SpatialDE Python package.
spatialDE(x, ...) ## S4 method for signature 'matrix' spatialDE(x, coordinates, verbose = FALSE) ## S4 method for signature 'SpatialExperiment' spatialDE(x, assay_type = "counts", verbose = FALSE)
spatialDE(x, ...) ## S4 method for signature 'matrix' spatialDE(x, coordinates, verbose = FALSE) ## S4 method for signature 'SpatialExperiment' spatialDE(x, assay_type = "counts", verbose = FALSE)
x |
A numeric Alternatively, a SpatialExperiment object. |
... |
For the generic, arguments to pass to specific methods. |
coordinates |
A For the SpatialExperiment method, coordinates are taken from
|
verbose |
A |
assay_type |
A |
A data.frame
with DE results where each row is a gene and columns
contain relevant statistics.
The most important columns are:
g
: the name of the gene
pval
: the p-value for spatial differential expression
qval
: the q-value, indicating significance after correcting for
multiple testing
l
: A parameter indicating the distance scale a gene changes expression
over
Davide Corso, Milan Malfait, Lambda Moses
Svensson, V., Teichmann, S. & Stegle, O. SpatialDE: identification of spatially variable genes. Nat Methods 15, 343–346 (2018). https://doi.org/10.1038/nmeth.4636
SpatialDE 1.1.3: the version of the Python package used under the hood.
The individual steps performed by this function: stabilize()
,
regress_out()
and run()
.
For further analysis of the DE results:
model_search()
and spatial_patterns()
.
## Mock up a SpatialExperiment object wit 400 cells and 3 genes set.seed(42) spe <- mockSVG(size = 20, tot_genes = 3, de_genes = 1, return_SPE = TRUE) ## Run spatialDE de_results <- spatialDE(spe) head(de_results)
## Mock up a SpatialExperiment object wit 400 cells and 3 genes set.seed(42) spe <- mockSVG(size = 20, tot_genes = 3, de_genes = 1, return_SPE = TRUE) ## Run spatialDE de_results <- spatialDE(spe) head(de_results)
Group spatially variable genes into spatial patterns using Automatic Expression Histology, using the SpatialDE Python package.
spatialPatterns(x, de_results, ...) ## S4 method for signature 'matrix' spatialPatterns( x, de_results, coordinates, qval_thresh = 0.05, n_patterns, length, verbose = FALSE ) ## S4 method for signature 'SpatialExperiment' spatialPatterns( x, de_results, qval_thresh = 0.05, n_patterns, length, assay_type = "counts", verbose = FALSE )
spatialPatterns(x, de_results, ...) ## S4 method for signature 'matrix' spatialPatterns( x, de_results, coordinates, qval_thresh = 0.05, n_patterns, length, verbose = FALSE ) ## S4 method for signature 'SpatialExperiment' spatialPatterns( x, de_results, qval_thresh = 0.05, n_patterns, length, assay_type = "counts", verbose = FALSE )
x |
A numeric Alternatively, a SpatialExperiment object. |
de_results |
|
... |
For the generic, arguments to pass to specific methods. |
coordinates |
A For the SpatialExperiment method, coordinates are taken from
|
qval_thresh |
|
n_patterns |
|
length |
|
verbose |
A |
assay_type |
A |
A list
of two data.frame
s (pattern_results, patterns):
pattern_results
: data.frame
with pattern membership information for each
gene.
patterns
the posterior mean underlying expression from genes in given
spatial patterns.
Davide Corso, Milan Malfait, Lambda Moses
Svensson, V., Teichmann, S. & Stegle, O. SpatialDE: identification of spatially variable genes. Nat Methods 15, 343–346 (2018). https://doi.org/10.1038/nmeth.4636
SpatialDE 1.1.3: the version of the Python package used under the hood.
The individual steps performed by this function: stabilize()
,
regress_out()
and spatial_patterns()
.
## Mock up a SpatialExperiment object wit 100 cells and 3 genes set.seed(42) spe <- mockSVG(size = 10, tot_genes = 3, de_genes = 1, return_SPE = TRUE) ## Run spatialDE de_results <- spatialDE(spe) spatial_patterns <- spatialPatterns(spe, de_results = de_results, qval_thresh = NULL, n_patterns = 4L, length = 1.5, verbose = FALSE ) head(spatial_patterns$pattern_results) head(spatial_patterns$patterns)
## Mock up a SpatialExperiment object wit 100 cells and 3 genes set.seed(42) spe <- mockSVG(size = 10, tot_genes = 3, de_genes = 1, return_SPE = TRUE) ## Run spatialDE de_results <- spatialDE(spe) spatial_patterns <- spatialPatterns(spe, de_results = de_results, qval_thresh = NULL, n_patterns = 4L, length = 1.5, verbose = FALSE ) head(spatial_patterns$pattern_results) head(spatial_patterns$patterns)
Stabilize variance of negative binomial data using Anscombe's approximation.
This function is a wrapper for stabilize
from the
NaiveDE Python package.
stabilize(counts)
stabilize(counts)
counts |
|
matrix
of variance stabilized counts.
## Mock up a SpatialExperiment object wit 400 cells and 3 genes set.seed(42) mock <- mockSVG(20, 3, 1) stabilized <- stabilize(mock$counts)
## Mock up a SpatialExperiment object wit 400 cells and 3 genes set.seed(42) mock <- mockSVG(20, 3, 1) stabilized <- stabilize(mock$counts)