Title: | Joint Dimension Reduction and Spatial Clustering |
---|---|
Description: | Joint dimension reduction and spatial clustering is conducted for Single-cell RNA sequencing and spatial transcriptomics data, and more details can be referred to Wei Liu, Xu Liao, Yi Yang, Huazhen Lin, Joe Yeong, Xiang Zhou, Xingjie Shi and Jin Liu. (2022) <doi:10.1093/nar/gkac219>. It is not only computationally efficient and scalable to the sample size increment, but also is capable of choosing the smoothness parameter and the number of clusters as well. |
Authors: | Wei Liu [aut, cre], Yi Yang [aut], Jin Liu [aut] |
Maintainer: | Wei Liu <[email protected]> |
License: | GPL-3 |
Version: | 3.4 |
Built: | 2024-11-16 03:21:10 UTC |
Source: | https://github.com/feiyoung/dr.sc |
A human dorsolateral prefrontal cortex dataset measured on the Visium platform, which invcludes 4634 spots and 500 genes, a subset of raw dataset at https://github.com/LieberInstitute/spatialLIBD.
nothing
Wei Liu
None
data("dlpfc151510")
data("dlpfc151510")
Joint dimension reduction and spatial clustering for scRNA-seq and spatial transcriptomics data
## S3 method for class 'Seurat' DR.SC(seu, K, q=15, platform= c('Visium', "ST", "Other_SRT", "scRNAseq"),...)
## S3 method for class 'Seurat' DR.SC(seu, K, q=15, platform= c('Visium', "ST", "Other_SRT", "scRNAseq"),...)
seu |
an object of class "Seurat". The details of this object are given under 'Details'. |
q |
a positive integer, specify the number of latent features to be extracted, default as 15. |
K |
a positive integer or integer vector, specify the number of clusters. When |
platform |
a string, specify the platform of the provided data, default as "Visium". There are more platforms to be chosen, including ("Visuim", "ST", "Other_SRT") and "scRNAseq", where the first group means there are spatial coordinates information in the metadata of seu, named "row" and "col" and a Hidden Markov random field is used to model the unobserved class label using spatial coordinates ("Other_SRT" represents the other SRT platforms except for 'Visium' and 'ST'), and the other group "scRNAseq" means there is no spatial information in object seu and a multinomial model is used to model the unobserved class labels. The platform helps to calculate the adjacency matrix. |
... |
Other arguments to pass into DR.SC_fit function. |
seu is an object named Seurat
, thich can easily created by R package Seurat.
DR-SC model can be applied to analyze both single cell RNA sequencing data and spatially resoved transcriptomics (SRT) data.
If the data are collected by the single cell RNA sequencing techonologies which means there is no spatial information in object seu then a multinomial model is used to model the unobserved class labels.
If the data is collected by the spatial transcriptomics technologies, then there are spatial coordinates information in the metadata of seu, named "row" and "col". DR-SC model uses a Hidden Markov random field to model the spatial coordinates. DR.SC
supports different platforms of SRT data, such as 'Visium', 'ST' and any other platforms 'Other_SRT'.
For lattice grids in ST platform (ST), the interior spot has four neighbors (left, right, up and down),the boundary spot has three neighbors, and the spot in the corner has two neighbors. For hexagon grids, such as spatial coordinate in 10X Visium platform (Visium), the interior spot has six neighbors. For the irregular coordinates in other platforms (Other_SRT), Euclidean distance is adopted to decide whether a spot is an neighbor of another spot. For example, if the Euclidean distance between spot A and B is less than a radius, then A is taken as the neighbourhood of B. See function getAdj for more details.
DR.SC returns a revised Seurat
object. There are two revisions in the seu: 1. the metadata is added a new column named spatial.drsc.cluster
that represents the clustering results from DR-SC model, and the Idents(seu)
is assigned with spatial.drsc.cluster
. 2. a DimReduc object named dr-sc
is added in the slot reductions
, which represents the features extracted by DR-SC model.
nothing
Wei Liu
None
## we generate the spatial transcriptomics data with lattice neighborhood, i.e. ST platform. seu <- gendata_RNAExp(height=10, width=10,p=50, K=4,platform="ST") library(Seurat) seu <- NormalizeData(seu, verbose=FALSE) # choose 100 highly variable features # seu <- FindVariableFeatures(seu, nfeatures = 100) # maxIter = 2 is only used for illustration, and user can use default. # seu1 <- DR.SC(seu, K=4, platform = 'ST', maxIter=2,verbose=FALSE) # choose spatially variable features (SVGs) seu <- FindSVGs(seu, nfeatures = 40, verbose=FALSE) # use SVGs to fit DR.SC model # maxIter = 2 is only used for illustration, and user can use default. seu1 <- DR.SC(seu, K=4,platform = 'ST', maxIter=2, verbose=TRUE)
## we generate the spatial transcriptomics data with lattice neighborhood, i.e. ST platform. seu <- gendata_RNAExp(height=10, width=10,p=50, K=4,platform="ST") library(Seurat) seu <- NormalizeData(seu, verbose=FALSE) # choose 100 highly variable features # seu <- FindVariableFeatures(seu, nfeatures = 100) # maxIter = 2 is only used for illustration, and user can use default. # seu1 <- DR.SC(seu, K=4, platform = 'ST', maxIter=2,verbose=FALSE) # choose spatially variable features (SVGs) seu <- FindSVGs(seu, nfeatures = 40, verbose=FALSE) # use SVGs to fit DR.SC model # maxIter = 2 is only used for illustration, and user can use default. seu1 <- DR.SC(seu, K=4,platform = 'ST', maxIter=2, verbose=TRUE)
Joint dimension reduction and spatial clustering for scRNA-seq and spatial transcriptomics data
DR.SC_fit(X, K, Adj_sp=NULL, q=15, error.heter= TRUE, beta_grid=seq(0.5, 5, by=0.5), maxIter=25, epsLogLik=1e-5, verbose=FALSE, maxIter_ICM=6, wpca.int=FALSE, int.model="EEE", approxPCA=FALSE, coreNum = 5)
DR.SC_fit(X, K, Adj_sp=NULL, q=15, error.heter= TRUE, beta_grid=seq(0.5, 5, by=0.5), maxIter=25, epsLogLik=1e-5, verbose=FALSE, maxIter_ICM=6, wpca.int=FALSE, int.model="EEE", approxPCA=FALSE, coreNum = 5)
X |
a sparse matrix with class |
K |
a positive integer allowing scalar or vector, specify the number of clusters in model fitting. |
Adj_sp |
an optional sparse matrix with class |
q |
a positive integer, specify the number of latent features to be extracted, default as 15. Usually, the choice of q is a trade-off between model complexity and fit to the data, and depends on the goals of the analysis and the structure of the data. A higher value will result in a more complex model with a higher number of parameters, which may lead to overfitting and poor generalization performance. On the other hand, a lower value will result in a simpler model with fewer parameters, but may also lead to underfitting and a poorer fit to the data. |
error.heter |
an optional logical value, whether use the heterogenous error for DR-SC model, default as |
beta_grid |
an optional vector of positive value, the candidate set of the smoothing parameter to be searched by the grid-search optimization approach. |
maxIter |
an optional positive value, represents the maximum iterations of EM. |
epsLogLik |
an optional positive vlaue, tolerance vlaue of relative variation rate of the observed pseudo log-loglikelihood value, defualt as '1e-5'. |
verbose |
an optional logical value, whether output the information of the ICM-EM algorithm. |
maxIter_ICM |
an optional positive value, represents the maximum iterations of ICM. |
wpca.int |
an optional logical value, means whether use the weighted PCA to obtain the initial values of loadings and other paramters, default as |
int.model |
an optional string, specify which Gaussian mixture model is used in evaluting the initial values for DR-SC, default as "EEE"; and see Mclust for more models' names. |
approxPCA |
an optional logical value, whether use approximated PCA to speed up the computation for initial values. |
coreNum |
an optional positive integer, means the number of thread used in parallel computating, default as 5. If the length of K is one, then coreNum will be set as 1 automatically. |
Nothing
DR.SC_fit returns a list with class "drscObject" with the following three components:
Objdrsc |
a list including the model fitting results, in which the number of elements is same as the length of K. |
out_param |
a numeric matrix used for model selection in MBIC. |
K_set |
a scalar or vector equal to input argument K. |
In addition, each element of "Objdrsc" is a list with the following comoponents:
cluster |
inferred class labels |
hZ |
extracted latent features. |
beta |
estimated smoothing parameter |
Mu |
mean vectors of mixtures components. |
Sigma |
covariance matrix of mixtures components. |
W |
estimated loading matrix |
Lam_vec |
estimated variance of errors in probabilistic PCA model |
loglik |
pseudo observed log-likelihood. |
nothing
Wei Liu
None
## we generate the spatial transcriptomics data with lattice neighborhood, i.e. ST platform. seu <- gendata_RNAExp(height=10, width=10,p=50, K=4) library(Seurat) seu <- NormalizeData(seu, verbose=FALSE) # choose 40 highly variable features using FindVariableFeatures in Seurat # seu <- FindVariableFeatures(seu, nfeatures = 40) # or choose 40 spatailly variable features using FindSVGs in DR.SC seu <- FindSVGs(seu, nfeatures = 40, verbose=FALSE) # users define the adjacency matrix Adj_sp <- getAdj(seu, platform = 'ST') if(class(seu@assays$RNA)=="Assay5"){ var.features <- seu@[email protected]$var.features var.features <- var.features[!is.na(var.features )] dat <- GetAssayData(seu, assay = "RNA", slot='data') X <- Matrix::t(dat[var.features,]) }else{ var.features <- seu@[email protected] X <- Matrix::t(seu[["RNA"]]@data[var.features,]) } # maxIter = 2 is only used for illustration, and user can use default. drscList <- DR.SC_fit(X,Adj_sp=Adj_sp, K=4, maxIter=2, verbose=TRUE)
## we generate the spatial transcriptomics data with lattice neighborhood, i.e. ST platform. seu <- gendata_RNAExp(height=10, width=10,p=50, K=4) library(Seurat) seu <- NormalizeData(seu, verbose=FALSE) # choose 40 highly variable features using FindVariableFeatures in Seurat # seu <- FindVariableFeatures(seu, nfeatures = 40) # or choose 40 spatailly variable features using FindSVGs in DR.SC seu <- FindSVGs(seu, nfeatures = 40, verbose=FALSE) # users define the adjacency matrix Adj_sp <- getAdj(seu, platform = 'ST') if(class(seu@assays$RNA)=="Assay5"){ var.features <- seu@assays$RNA@meta.data$var.features var.features <- var.features[!is.na(var.features )] dat <- GetAssayData(seu, assay = "RNA", slot='data') X <- Matrix::t(dat[var.features,]) }else{ var.features <- seu@assays$RNA@var.features X <- Matrix::t(seu[["RNA"]]@data[var.features,]) } # maxIter = 2 is only used for illustration, and user can use default. drscList <- DR.SC_fit(X,Adj_sp=Adj_sp, K=4, maxIter=2, verbose=TRUE)
Intuitive way of visualizing how cell types changes across the embeddings obatined by DR-SC.
drscPlot(seu, dims=1:5, visu.method='tSNE',...)
drscPlot(seu, dims=1:5, visu.method='tSNE',...)
seu |
an object of class "Seurat" obtained by DR.SC. |
dims |
a positive integer to specify the number of latent features for visualiztion. |
visu.method |
a string including 'tSNE' or "UMAP". |
... |
Other arguments passing to DimPlot function. |
Nothing
return a ggplot2 object.
nothing
Wei Liu
None
None
## we generate the spatial transcriptomics data with lattice neighborhood, i.e. ST platform. seu <- gendata_RNAExp(height=10, width=10,p=50, K=4) library(Seurat) seu <- NormalizeData(seu) # choose spatially variable features seu <- FindSVGs(seu) # use SVGs to fit DR.SC model # maxIter = 2 is only used for illustration, and user can use default. seu1 <- DR.SC(seu, K=4,platform = 'ST', maxIter = 2,verbose=FALSE) drscPlot(seu1)
## we generate the spatial transcriptomics data with lattice neighborhood, i.e. ST platform. seu <- gendata_RNAExp(height=10, width=10,p=50, K=4) library(Seurat) seu <- NormalizeData(seu) # choose spatially variable features seu <- FindSVGs(seu) # use SVGs to fit DR.SC model # maxIter = 2 is only used for illustration, and user can use default. seu1 <- DR.SC(seu, K=4,platform = 'ST', maxIter = 2,verbose=FALSE) drscPlot(seu1)
Identifies features that have spatially variation along spots using SPARK-X.
FindSVGs(seu, nfeatures=2000, covariates=NULL, preHVGs=5000,num_core=1, verbose=TRUE)
FindSVGs(seu, nfeatures=2000, covariates=NULL, preHVGs=5000,num_core=1, verbose=TRUE)
seu |
an object of class "Seurat". |
nfeatures |
a positive integer, means how many spatially variable genes to be chosen. If there are less than 2000 features in seu, then all features are identified. |
covariates |
a covariate matrix named control variable matrix whose number of rows is equal to the number of columns of seu. |
preHVGs |
a positive integer, the number of highly variable genes selected for speeding up computation of SPARK-X in selecting spatially variable features. |
num_core |
an optional positive integer, specify the cores used for identifying the SVGs in parallel. |
verbose |
an optional logical value, whether output the related information. |
Nothing
return a revised Seurat object by adding three columns named "is.SVGs", "order.SVGs" and "adjusted.pval.SVGs" in the meta.features of default Assay.
nothing
Zhu, J., Sun, S., Zhou, X.: Spark-x: non-parametric modeling enables scalable and robust detection of spatialexpression patterns for large spatial transcriptomic studies. Genome Biology 22(1), 1-25 (2021)
seu<-gendata_RNAExp(height=20, width=20,p=200, K=4) seu<-FindSVGs(seu, nfeatures=100) topSVGs(seu)
seu<-gendata_RNAExp(height=20, width=20,p=200, K=4) seu<-FindSVGs(seu, nfeatures=100) topSVGs(seu)
Generate simulated spatial transcriptomics data or scRNAseq data.
gendata_RNAExp(height=30, width=30, platform="ST", p =100, q=10, K=7, G=4,sigma2=1, tau=8, seed=1, view=FALSE)
gendata_RNAExp(height=30, width=30, platform="ST", p =100, q=10, K=7, G=4,sigma2=1, tau=8, seed=1, view=FALSE)
height , width
|
Height and width of lattice grids for generating spatial coordinates. n=height * width cells for scRNAseq data |
platform |
set the platform for the simulated data, only support 'ST' and 'scRNAseq'. |
p |
number of genes to generate. |
q |
number of true latent features to generate gene expression |
K |
number of clusters (cell types). |
seed |
random seed for generate data |
G |
the number of neighbors. The latter must be one of G = 4 or G = 8, which respectively correspond to a first order and a second order dependency structure. By default, G = 4. |
sigma2 |
Variance of error term in probabilitic PCA model. |
tau |
a positive factor of mixture mean values. |
view |
Logical value indicating whether the draw should be printed. Do not display the optional borders. |
Nothing
return a "Seurat" object. If platform="ST"
, then the metadata of this Seurat object will include two columns with names "row" and "col" which are the spatial coordinates; If platform="scRNAseq"
, then the metadata of this Seurat object will not have them.
nothing
Wei Liu
None
None
## we generate the spatial transcriptomics data with lattice neighborhood, i.e. ST platform. seu <- gendata_RNAExp(height=20, width=20,p=200, K=4) seu ## generate scRNAseq data seu <- gendata_RNAExp(height=20, width=20, platform="scRNAseq", p=100, K=4) seu
## we generate the spatial transcriptomics data with lattice neighborhood, i.e. ST platform. seu <- gendata_RNAExp(height=20, width=20,p=200, K=4) seu ## generate scRNAseq data seu <- gendata_RNAExp(height=20, width=20, platform="scRNAseq", p=100, K=4) seu
Calculate the adjacency matrix for the spatial transcriptomics data measured on 10X Visium, ST or other platforms as a Seurat object.
## S3 method for class 'Seurat' getAdj(obj, platform = c('Visium', "ST", "Other_SRT"), ...)
## S3 method for class 'Seurat' getAdj(obj, platform = c('Visium', "ST", "Other_SRT"), ...)
obj |
an object with class "Seurat", there are spatial coordinates information in the metadata of obj, named "row" and "col", where first column is x-axis coordinate, the second column is y-axis coordinate. getAdj_manual and getAdj_auto supports multi-dimensional spatial coordinates with a matrix as input. |
platform |
a string, specify the platform of the provided data, default as "Visium". There are more platforms to be chosen, including ("Visuim", "ST", "Other_SRT"), which means there are spatial coordinates information in the metadata of obj, named "row" and "col". The platform helps to calculate the adjacency matrix by defining the neighborhoods. |
... |
Other arguments to pass into getAdj_auto function. |
For lattice grids in ST platform (ST), the interior spot has four neighbors (left, right, up and down),the boundary spot has three neighbors, and the spot in the corner has two neighbors. For hexagon grids, such as spatial coordinate in 10X Visium platform (Visium), the interior spot has six neighbors. For the irregular coordinates in other platforms (Other_SRT), Euclidean distance is adopted to decide whether a spot is an neighbor of another spot. For example, if the Euclidean distance between spot A and B is less than a radius, then A is taken as the neighbourhood of B. See functions getAdj_auto and getAdj_manual for more details.
Return a dgCMatrix
object recording the information of neighborhoods about each spot.
nothing
Wei Liu
## S3 method for class "Seurat" seu <- gendata_RNAExp(height=20, width=20,p=200, K=4) Adj_sp <- getAdj(seu, platform = 'ST')
## S3 method for class "Seurat" seu <- gendata_RNAExp(height=20, width=20,p=200, K=4) Adj_sp <- getAdj(seu, platform = 'ST')
an efficient function to find the radius by bi-section method , then find neighbors based on the matrix of position, which ensures that each spot has approximately lower.med~upper.med neighbors in the sense of median.
getAdj_auto(pos, lower.med=4, upper.med=6, radius.upper= NULL)
getAdj_auto(pos, lower.med=4, upper.med=6, radius.upper= NULL)
pos |
a n-by-2 matrix of position. |
lower.med |
an integer, the lower bound of median number of neighbors among all spots. |
upper.med |
an integer, the upper bound of median number of neighbors among all spots. |
radius.upper |
a real, the upper bound of radius, default as NULL. If |
A sparse adjacency matrix containing the neighbourhood.
an efficient function to find the neighbors based on the matrix of position and a pre-defined radius.
getAdj_manual(pos, radius)
getAdj_manual(pos, radius)
pos |
is a n-by-d matrix of position, where n is the number of spots, and d is the dimension of coordinates. |
radius |
is a threashold of Euclidean distance to decide whether a spot is an neighborhood of another spot. For example, if the Euclidean distance between spot A and B is less than radius, then A is taken as the neighbourhood of B. |
A sparse matrix containing the neighbourhood
an efficient function to find the neighborhood based on the matrix of position and a pre-defined cutoff
getneighborhood_fast(x, radius)
getneighborhood_fast(x, radius)
x |
is a n-by-2 matrix of position. |
radius |
is a threashold of Euclidean distance to decide whether a spot is an neighborhood of another spot. For example, if the Euclidean distance between spot A and B is less than cutoff, then A is taken as the neighbourhood of B. |
A sparse matrix containing the neighbourhood
Intuitive way of visualizing how modified BIC values changes across different number of clusters
mbicPlot(seu, criteria="MBIC")
mbicPlot(seu, criteria="MBIC")
seu |
an object of class |
criteria |
a string specifying the information criteria such as AIC, BIC and MBIC to be plotted, default as MBIC. |
Nothing
return a ggplot2 object.
nothing
Wei Liu
None
None
## we generate the spatial transcriptomics data with lattice neighborhood, i.e. ST platform. seu <- gendata_RNAExp(height=20, width=20,p=100, K=4) library(Seurat) seu <- NormalizeData(seu) # choose spatially variable features seu <- FindSVGs(seu) ## Just for illustrating the usage of mbicPlot seu[["RNA"]]@misc[['icMat']] <- data.frame(K=2:5, MBIC=c(105, 101, 99, 108)) mbicPlot(seu)
## we generate the spatial transcriptomics data with lattice neighborhood, i.e. ST platform. seu <- gendata_RNAExp(height=20, width=20,p=100, K=4) library(Seurat) seu <- NormalizeData(seu) # choose spatially variable features seu <- FindSVGs(seu) ## Just for illustrating the usage of mbicPlot seu[["RNA"]]@misc[['icMat']] <- data.frame(K=2:5, MBIC=c(105, 101, 99, 108)) mbicPlot(seu)
Read the spatial transcriptomics data measured on 10X Visium platform as a Seurat object, where the spatial coordinates are saved in the metadata, named "row" and "col".
read10XVisium(dirname)
read10XVisium(dirname)
dirname |
A string, the dictory of Visium datasets |
Nothing
return a Seurat object.
nothing
Wei Liu
None
None
## Not run: ## set your file directory, then read it. data_name <- "D/HCC" HCC1 <- read10XVisium(data_name) ## End(Not run)
## Not run: ## set your file directory, then read it. data_name <- "D/HCC" HCC1 <- read10XVisium(data_name) ## End(Not run)
Read the single cell RNA sequencing data measured on scRNA sequencing platform as a Seurat object.
readscRNAseq(mtx, cells, features, ...)
readscRNAseq(mtx, cells, features, ...)
mtx |
a string, ame or remote URL of the mtx file |
cells |
a string, Name or remote URL of the cells/barcodes file |
features |
a string, Name or remote URL of the features/genes file |
... |
the arguments passing to ReadMtx |
Nothing
return a Seurat object including expression matrix.
nothing
Wei Liu
None
None
## Not run: ### set the file directory, then read it. seu <- readscRNAseq(mtx="GSM3755564_16_Liver_Treg_matrix.mtx.gz", features='GSM3755564_16_Liver_Treg_genes.tsv.gz', cells='GSM3755564_16_Liver_Treg_barcodes.tsv.gz' ) seu ## End(Not run)
## Not run: ### set the file directory, then read it. seu <- readscRNAseq(mtx="GSM3755564_16_Liver_Treg_matrix.mtx.gz", features='GSM3755564_16_Liver_Treg_genes.tsv.gz', cells='GSM3755564_16_Liver_Treg_barcodes.tsv.gz' ) seu ## End(Not run)
Run a weighted PCA dimensionality reduction
RunWPCA(object, q=15) ### S3 method for class "Seurat" ## RunWPCA(object, q=15) ### S3 method for class "matrix" ## RunWPCA(object, q=15) ### S3 method for class "dgCMatrix" ## RunWPCA(object, q=15)
RunWPCA(object, q=15) ### S3 method for class "Seurat" ## RunWPCA(object, q=15) ### S3 method for class "matrix" ## RunWPCA(object, q=15) ### S3 method for class "dgCMatrix" ## RunWPCA(object, q=15)
object |
an object named "Seurat", "maxtrix" or "dgCMatrix". The object of class "Seurat" must include slot "scale.data". |
q |
an optional positive integer, specify the number of features to be extracted. |
Nothing
For Seurat object, return a Seurat object. For objcet "matrix" and "dgCMatrix", return a object "matrix" with rownames same as the colnames of X
, and colnames "WPCA1" to "WPCAq".
nothing
Wei Liu
Bai, J. and Liao, Y. (2017). Inferences in panel data with interactive effects using large covariance matrices. Journal of Econometrics, 200(1):59–78.
None
## Not run: library(Seurat) seu <- gendata_RNAExp(height=20, width=20,p=100, K=4) ## log-normalization seu <- NormalizeData(seu) ## seu <- FindVariableFeatures(seu, nfeatures=80) ## Scale seu <- ScaleData(seu) ## Run WPCA seu <- RunWPCA(seu) seu ## Run tSNE based on wpca seu <- RunTSNE(seu, reduction='wpca') seu ## Find SVGs seu <- FindSVGs(seu, nfeatures=80) (genes <- topSVGs(seu, ntop=10)) Idents(seu) <- factor(paste0("cluster", seu$true_clusters), levels=paste0("cluster",1:4)) RidgePlot(seu, features = genes[1:2], ncol = 2) FeaturePlot(seu, features = genes[1:2], reduction = 'tsne' ,ncol=2) ## End(Not run)
## Not run: library(Seurat) seu <- gendata_RNAExp(height=20, width=20,p=100, K=4) ## log-normalization seu <- NormalizeData(seu) ## seu <- FindVariableFeatures(seu, nfeatures=80) ## Scale seu <- ScaleData(seu) ## Run WPCA seu <- RunWPCA(seu) seu ## Run tSNE based on wpca seu <- RunTSNE(seu, reduction='wpca') seu ## Find SVGs seu <- FindSVGs(seu, nfeatures=80) (genes <- topSVGs(seu, ntop=10)) Idents(seu) <- factor(paste0("cluster", seu$true_clusters), levels=paste0("cluster",1:4)) RidgePlot(seu, features = genes[1:2], ncol = 2) FeaturePlot(seu, features = genes[1:2], reduction = 'tsne' ,ncol=2) ## End(Not run)
Select the number of clusters by specified criteria.
selectModel(obj, criteria = 'MBIC', pen.const=1) ## S3 method for class 'drscObject' selectModel(obj, criteria = 'MBIC', pen.const=1) ## S3 method for class 'Seurat' selectModel(obj, criteria = 'MBIC', pen.const=1)
selectModel(obj, criteria = 'MBIC', pen.const=1) ## S3 method for class 'drscObject' selectModel(obj, criteria = 'MBIC', pen.const=1) ## S3 method for class 'Seurat' selectModel(obj, criteria = 'MBIC', pen.const=1)
S
obj |
an object with class |
criteria |
a string, specify the criteria used for selecting the number of clusters, supporting "MBIC", "BIC" and "AIC". |
pen.const |
an optional positive value, the adjusted constant used in the MBIC criteria. It usually takes value between 0.1 to 1. |
For S3 method of Seurat
, it return a revised "Seurat" object with updated Idents(seu)
, spatial.drsc.cluster
in the metadata and DimReduc object named dr-sc
in the slot reductions
. For S3 method of drscObject
, it returns a list with the following components:
bestK |
the selected number of clusters. |
cluster |
inferred class labels |
hZ |
extracted latent features. |
icMat |
a numeric matrix including the criteria value for each number of clusters K. |
nothing
Wei Liu
seu <- gendata_RNAExp(height=10, width=10,p=50, K=4) library(Seurat) seu <- NormalizeData(seu, verbose=FALSE) # or choose 40 spatailly variable features using FindSVGs in DR.SC seu <- FindSVGs(seu, nfeatures = 40, verbose=FALSE) # users define the adjacency matrix Adj_sp <- getAdj(seu, platform = 'ST') dat <- GetAssayData(seu, assay = "RNA", slot='data') X <- Matrix::t(dat) # maxIter = 2 is only used for illustration, and user can use default. drscList <- DR.SC_fit(X,Adj_sp=Adj_sp ,K=4, maxIter=2, verbose=TRUE) drsc1 <- selectModel(drscList) str(drsc1)
seu <- gendata_RNAExp(height=10, width=10,p=50, K=4) library(Seurat) seu <- NormalizeData(seu, verbose=FALSE) # or choose 40 spatailly variable features using FindSVGs in DR.SC seu <- FindSVGs(seu, nfeatures = 40, verbose=FALSE) # users define the adjacency matrix Adj_sp <- getAdj(seu, platform = 'ST') dat <- GetAssayData(seu, assay = "RNA", slot='data') X <- Matrix::t(dat) # maxIter = 2 is only used for illustration, and user can use default. drscList <- DR.SC_fit(X,Adj_sp=Adj_sp ,K=4, maxIter=2, verbose=TRUE) drsc1 <- selectModel(drscList) str(drsc1)
Calculate column-wise or row-wise mean
sp_means_Rcpp(sp_data, rowMeans = FALSE)
sp_means_Rcpp(sp_data, rowMeans = FALSE)
sp_data |
A sparse matrix |
rowMeans |
A boolean value, whether to calculate row-wise mean |
A n x 1 or p x 1 matrix
Calculate column-wise or row-wise sum
sp_sums_Rcpp(sp_data, rowSums = FALSE)
sp_sums_Rcpp(sp_data, rowSums = FALSE)
sp_data |
A sparse matrix |
rowSums |
A boolean value, whether to calculate row-wise sum |
A n x 1 or p x 1 matrix
Intuitive way of visualizing how cell types changes across the spatial locations.
spatialPlotClusters(seu)
spatialPlotClusters(seu)
seu |
an object of class "Seurat" obtained by DR.SC. |
Nothing
return a ggplot2 object.
nothing
Wei Liu
None
None
## we generate the spatial transcriptomics data with lattice neighborhood, i.e. ST platform. seu <- gendata_RNAExp(height=10, width=10,p=50, K=4) library(Seurat) seu <- NormalizeData(seu) # choose spatially variable features using Seurat seu <- FindSVGs(seu) # use SVGs to fit DR.SC model # maxIter = 2 is only used for illustration, and user can use default. seu1 <- DR.SC(seu, K=4,platform = 'ST', maxIter=2,verbose=FALSE) spatialPlotClusters(seu1)
## we generate the spatial transcriptomics data with lattice neighborhood, i.e. ST platform. seu <- gendata_RNAExp(height=10, width=10,p=50, K=4) library(Seurat) seu <- NormalizeData(seu) # choose spatially variable features using Seurat seu <- FindSVGs(seu) # use SVGs to fit DR.SC model # maxIter = 2 is only used for illustration, and user can use default. seu1 <- DR.SC(seu, K=4,platform = 'ST', maxIter=2,verbose=FALSE) spatialPlotClusters(seu1)
Return top n spatially variable genes given a Seurat object performed by FindSVGs.
topSVGs(seu, ntop=5)
topSVGs(seu, ntop=5)
seu |
an object of class "Seurat". |
ntop |
an optional positive integer, means how many spatially variable genes to access. |
Nothing
return a character vector including the names of SVGs.
nothing
Wei Liu
None
seu <- gendata_RNAExp(height=20, width=20,p=200, K=4) seu <- FindSVGs(seu, nfeatures=100, verbose=FALSE) (genes <- topSVGs(seu, ntop=10))
seu <- gendata_RNAExp(height=20, width=20,p=200, K=4) seu <- FindSVGs(seu, nfeatures=100, verbose=FALSE) (genes <- topSVGs(seu, ntop=10))