Title: | Thin plate models to detect spatially variable genes |
---|---|
Description: | The goal of `tpSVG` is to detect and visualize spatial variation in the gene expression for spatially resolved transcriptomics data analysis. Specifically, `tpSVG` introduces a family of count-based models, with generalizable parametric assumptions such as Poisson distribution or negative binomial distribution. In addition, comparing to currently available count-based model for spatially resolved data analysis, the `tpSVG` models improves computational time, and hence greatly improves the applicability of count-based models in SRT data analysis. |
Authors: | Boyi Guo [aut, cre] , Lukas M. Weber [ctb] , Stephanie C. Hicks [aut] |
Maintainer: | Boyi Guo <[email protected]> |
License: | MIT + file LICENSE |
Version: | 1.3.0 |
Built: | 2024-11-19 04:44:24 UTC |
Source: | https://github.com/bioc/tpSVG |
Thin Plate Spline Model to Detect Spatially Variable Genes
tpSVG( input, spatial_coords = NULL, X = NULL, family = poisson(), offset = log(input$sizeFactor), weights = NULL, assay_name = "counts", n_threads = 1, BPPARAM = NULL, verbose = FALSE, ... )
tpSVG( input, spatial_coords = NULL, X = NULL, family = poisson(), offset = log(input$sizeFactor), weights = NULL, assay_name = "counts", n_threads = 1, BPPARAM = NULL, verbose = FALSE, ... )
input |
|
spatial_coords |
|
X |
|
family |
a description of the error distribution and link function
to be used in the model. Currently support two distributions |
offset |
This can be used to account for technician variation when
|
weights |
Reserved for future development, e.g. correcting mean-var relationship for Gaussian models. Please use with caution. |
assay_name |
|
n_threads |
|
BPPARAM |
|
verbose |
|
... |
Reserved for future arguments. |
If the input was provided as a SpatialExperiment
object, the
output values are returned as additional columns in the rowData
slot
of the input object. If the input was provided as a numeric
matrix
of values, the output is returned as a numeric
matrix. The output
values include p-values without any adjustment and statistics reporting
reporting the thinplate spline model. The test_stat
entry of the
returned object is the test statistic for the corresponding model,
that is F statistics for the gaussian model and the Chi-squared statistics
for generalized models.
library(SpatialExperiment) library(STexampleData) library(scran) library(nnSVG) # load example dataset from STexampleData package spe <- Visium_humanDLPFC() # preprocessing steps # keep only spots over tissue spe <- spe[, colData(spe)$in_tissue == 1] # skip spot-level quality control, since this has been performed previously # on this dataset # Add library size spe <- addPerCellQCMetrics(spe) # filter low-expressed and mitochondrial genes spe <- filter_genes(spe) # calculate logcounts (log-transformed normalized counts) using scran package # using library size factors spe <- computeLibraryFactors(spe) spe <- logNormCounts(spe) # select small number of genes for faster runtime in this example set.seed(123) ix <- sample(seq_len(nrow(spe)), 4) spe <- spe[ix, ] # run tpSVG set.seed(123) # Gaussian Model spe_gaus <- tpSVG( spe, family = gaussian(), assay_name = "logcounts" ) # Poisson Model spe_poisson <- tpSVG( spe, family = poisson, assay_name = "counts", offset = log(spe$sizeFactor) # Natural log library size )
library(SpatialExperiment) library(STexampleData) library(scran) library(nnSVG) # load example dataset from STexampleData package spe <- Visium_humanDLPFC() # preprocessing steps # keep only spots over tissue spe <- spe[, colData(spe)$in_tissue == 1] # skip spot-level quality control, since this has been performed previously # on this dataset # Add library size spe <- addPerCellQCMetrics(spe) # filter low-expressed and mitochondrial genes spe <- filter_genes(spe) # calculate logcounts (log-transformed normalized counts) using scran package # using library size factors spe <- computeLibraryFactors(spe) spe <- logNormCounts(spe) # select small number of genes for faster runtime in this example set.seed(123) ix <- sample(seq_len(nrow(spe)), 4) spe <- spe[ix, ] # run tpSVG set.seed(123) # Gaussian Model spe_gaus <- tpSVG( spe, family = gaussian(), assay_name = "logcounts" ) # Poisson Model spe_poisson <- tpSVG( spe, family = poisson, assay_name = "counts", offset = log(spe$sizeFactor) # Natural log library size )