Title: | stJoincount - Join count statistic for quantifying spatial correlation between clusters |
---|---|
Description: | stJoincount facilitates the application of join count analysis to spatial transcriptomic data generated from the 10x Genomics Visium platform. This tool first converts a labeled spatial tissue map into a raster object, in which each spatial feature is represented by a pixel coded by label assignment. This process includes automatic calculation of optimal raster resolution and extent for the sample. A neighbors list is then created from the rasterized sample, in which adjacent and diagonal neighbors for each pixel are identified. After adding binary spatial weights to the neighbors list, a multi-categorical join count analysis is performed to tabulate "joins" between all possible combinations of label pairs. The function returns the observed join counts, the expected count under conditions of spatial randomness, and the variance calculated under non-free sampling. The z-score is then calculated as the difference between observed and expected counts, divided by the square root of the variance. |
Authors: | Jiarong Song [cre, aut] , Rania Bassiouni [aut], David Craig [aut] |
Maintainer: | Jiarong Song <[email protected]> |
License: | MIT + file LICENSE |
Version: | 1.9.0 |
Built: | 2024-11-19 04:29:22 UTC |
Source: | https://github.com/bioc/stJoincount |
Create a dictionary with categorical cluster labels as values and their converted numerical labels as keys
customDict(sampleInfo)
customDict(sampleInfo)
sampleInfo |
A data.frame contains the pixel information and cluster labels for each barcode of a human breast cancer sample. The index contains barcodes, and at least three other columns that have these information are required and the column names should be the same as following: "imagerow": The row pixel coordinate of the center of the spot "imagecol": The column pixel coordinate of the center of the spot "Cluster": The label that corresponding to this barcode |
A dictionary with categorical cluster labels as values and their converted numerical labels as keys.
fpath <- system.file("extdata", "dataframe.rda", package="stJoincount") load(fpath) nameList <- customDict(humanBC)
fpath <- system.file("extdata", "dataframe.rda", package="stJoincount") load(fpath) nameList <- customDict(humanBC)
This data preparation function creates a data.frame from a Seurat Object. It extracts the pixel information and cluster labels of each barcode from user's input Seurat Object and generate a data.frame with a certain format which is required for the algorithm. If the user has customized labels, this function will change the column name to "Cluster" when generating the data.frame to make it consistent to the required format.
dataPrepFromSeurat(SeuratObj, label)
dataPrepFromSeurat(SeuratObj, label)
SeuratObj |
input Seurat object that contains labels for each barcode |
label |
the column name of the label information in "meta.data" |
A data.frame contains the pixel information and cluster labels for each barcode the sample. The index contains barcodes, and at least three other columns that have these information are required and the column names should be the same as following: "imagerow": The row pixel coordinate of the center of the spot "imagecol": The column pixel coordinate of the center of the spot "Cluster": The label that corresponding to this barcode
fpath <- system.file("extdata", "SeuratBC.rda", package="stJoincount") load(fpath) df <- dataPrepFromSeurat(seuratBC, "Cluster")
fpath <- system.file("extdata", "SeuratBC.rda", package="stJoincount") load(fpath) df <- dataPrepFromSeurat(seuratBC, "Cluster")
This data preparation function creates a data.frame form Spatial Experiment Object. It extracts the pixel information and cluster labels of each barcode from user's input Seurat Object and generate a data.frame with a certain format which is required for the algorithm. If the user has customized labels, this function will change the column name to "Cluster" when generating the data.frame to make it consist to the required format.
dataPrepFromSpE(SpeObj, label)
dataPrepFromSpE(SpeObj, label)
SpeObj |
input SpatialExperiment object that contains labels for each barcode |
label |
the column name of the label information in "colData" |
A data.frame contains the pixel information and cluster labels for each barcode of the sample. The index contains barcodes, and at least three other columns that have these information are required and the column names should be the same as following: "imagerow": The row pixel coordinate of the center of the spot "imagecol": The column pixel coordinate of the center of the spot "Cluster": The label that corresponding to this barcode
fpath <- system.file("extdata", "SpeBC.rda", package="stJoincount") load(fpath) df <- dataPrepFromSpE(SpeObjBC, "label")
fpath <- system.file("extdata", "SpeBC.rda", package="stJoincount") load(fpath) df <- dataPrepFromSpE(SpeObjBC, "label")
When we create the rasterlayer, there will be a rectangular range. It is often necessary to provide a buffer to ensure that subsequent functions do not result in blank or missed pixels. This function is to find the right buffer for the sample coordinates so that each cluster is not lost in the process of converting a spot to a pixel.
extentBuffer(sampleInfo)
extentBuffer(sampleInfo)
sampleInfo |
A data.frame contains the pixel information and cluster labels for each barcode of a human breast cancer sample. The index contains barcodes, and at least three other columns that have these information are required and the column names should be the same as following: "imagerow": The row pixel coordinate of the center of the spot "imagecol": The column pixel coordinate of the center of the spot "Cluster": The label that corresponding to this barcode This data.frame can be produced by "dataPrepFromSeurat()/dataPrepFromSpE" functions |
optimal number of buffer for extent
fpath <- system.file("extdata", "dataframe.rda", package="stJoincount") load(fpath) n <- extentBuffer(humanBC)
fpath <- system.file("extdata", "dataframe.rda", package="stJoincount") load(fpath) n <- extentBuffer(humanBC)
This function performes multi-categorical join count analysis of the rasterized sample. A neighbors list is then created from the rasterized sample, in which adjacent and diagonal neighbors for each pixel are identified.
joincountAnalysis(mosaicIntegration)
joincountAnalysis(mosaicIntegration)
mosaicIntegration |
A raster object converted from a labeled spatial tissue map from Function rasterization. |
A data.frame that contains the observed join counts, the expected count under conditions of spatial randomness, the variance calculated under non-free sampling, and calculated Z-score.
fpath <- system.file("extdata", "dataframe.rda", package="stJoincount") load(fpath) mosaicIntegration <- rasterizeEachCluster(humanBC) joincount.result <- joincountAnalysis(mosaicIntegration)
fpath <- system.file("extdata", "dataframe.rda", package="stJoincount") load(fpath) mosaicIntegration <- rasterizeEachCluster(humanBC) joincount.result <- joincountAnalysis(mosaicIntegration)
Visualization of the rasterization results and label coding of the sample.
mosaicIntPlot(sampleInfo, mosaicIntegration)
mosaicIntPlot(sampleInfo, mosaicIntegration)
sampleInfo |
A dataset of a human breast cancer sample containing the pixel information and cluster labels for each barcode. |
mosaicIntegration |
A raster object converted from a labeled spatial tissue map. |
A mosaic plot with labeled pixels.
fpath <- system.file("extdata", "dataframe.rda", package="stJoincount") load(fpath) mosaicIntegration <- rasterizeEachCluster(humanBC) mosaicIntPlot(humanBC, mosaicIntegration)
fpath <- system.file("extdata", "dataframe.rda", package="stJoincount") load(fpath) mosaicIntegration <- rasterizeEachCluster(humanBC) mosaicIntPlot(humanBC, mosaicIntegration)
Converts a labeled spatial tissue map into a raster object, in which each spatial cluster is represented by a pixel coded by label assignment.
rasterizeEachCluster(sampleInfo)
rasterizeEachCluster(sampleInfo)
sampleInfo |
A data.frame contains the pixel information and cluster labels for each barcode of a human breast cancer sample. The index contains barcodes, and at least three other columns that have these information are required and the column names should be the same as following: "imagerow": The row pixel coordinate of the center of the spot "imagecol": The column pixel coordinate of the center of the spot "Cluster": The label that corresponding to this barcode |
This function returns a class of RasterLayer. This raster object is converted from a labeled spatial tissue map.
fpath <- system.file("extdata", "dataframe.rda", package="stJoincount") load(fpath) mosaicIntegration <- rasterizeEachCluster(humanBC)
fpath <- system.file("extdata", "dataframe.rda", package="stJoincount") load(fpath) mosaicIntegration <- rasterizeEachCluster(humanBC)
When sample coordinates finds a suitable buffer to ensure that each cluster is not lost in the process of converting the spot to pixel, apply this buffer to this function to find a suitable rectangle for the rasterlayer
rasterPrep(sampleInfo, n)
rasterPrep(sampleInfo, n)
sampleInfo |
A data.frame contains the pixel information and cluster labels for each barcode of a human breast cancer sample. The index contains barcodes, and at least three other columns that have these information are required and the column names should be the same as following: "imagerow": The row pixel coordinate of the center of the spot "imagecol": The column pixel coordinate of the center of the spot "Cluster": The label that corresponding to this barcode |
n |
buffer for extent (from function extentBuffer). |
This function returns a class of RasterLayer. This is a raster layer with calculated resolution and extent with buffer applied
fpath <- system.file("extdata", "dataframe.rda", package="stJoincount") load(fpath) raster <- rasterPrep(humanBC, 15)
fpath <- system.file("extdata", "dataframe.rda", package="stJoincount") load(fpath) raster <- rasterPrep(humanBC, 15)
Automatic calculation of optimal raster resolution for the sample.
resolutionCalc(sampleInfo)
resolutionCalc(sampleInfo)
sampleInfo |
A data.frame contains the pixel information and cluster labels for each barcode of a human breast cancer sample. The index contains barcodes, and at least three other columns that have these information are required and the column names should be the same as following: "imagerow": The row pixel coordinate of the center of the spot "imagecol": The column pixel coordinate of the center of the spot "Cluster": The label that corresponding to this barcode |
A list that contains length and height of resolution.
fpath <- system.file("extdata", "dataframe.rda", package="stJoincount") load(fpath) resolutionList <- resolutionCalc(humanBC)
fpath <- system.file("extdata", "dataframe.rda", package="stJoincount") load(fpath) resolutionList <- resolutionCalc(humanBC)
This function provides a heatmap of z-scores resulting from the join count analysis for all possible label pairs.
zscoreMatrix(sampleInfo, joincount.result)
zscoreMatrix(sampleInfo, joincount.result)
sampleInfo |
A data.frame contains the pixel information and cluster labels for each barcode of a human breast cancer sample. The index contains barcodes, and at least three other columns that have these information are required and the column names should be the same as following: "imagerow": The row pixel coordinate of the center of the spot "imagecol": The column pixel coordinate of the center of the spot "Cluster": The label that corresponding to this barcode |
joincount.result |
calculated result from join count analysis |
A data.frame that has a z-score matrix resulting from the join count analysis for all possible label pairs
fpath <- system.file("extdata", "dataframe.rda", package="stJoincount") load(fpath) mosaicIntegration <- rasterizeEachCluster(humanBC) joincount.result <- joincountAnalysis(mosaicIntegration) matrix <- zscoreMatrix(humanBC, joincount.result)
fpath <- system.file("extdata", "dataframe.rda", package="stJoincount") load(fpath) mosaicIntegration <- rasterizeEachCluster(humanBC) joincount.result <- joincountAnalysis(mosaicIntegration) matrix <- zscoreMatrix(humanBC, joincount.result)
Visulization of Z-score heatmap.
zscorePlot(zscoreMatrix)
zscorePlot(zscoreMatrix)
zscoreMatrix |
calculated and reshaped z-score matirx from join count analysis. |
A Heatmap plot
fpath <- system.file("extdata", "dataframe.rda", package="stJoincount") load(fpath) mosaicIntegration <- rasterizeEachCluster(humanBC) joincount.result <- joincountAnalysis(mosaicIntegration) matrix <- zscoreMatrix(humanBC, joincount.result) zscorePlot(matrix)
fpath <- system.file("extdata", "dataframe.rda", package="stJoincount") load(fpath) mosaicIntegration <- rasterizeEachCluster(humanBC) joincount.result <- joincountAnalysis(mosaicIntegration) matrix <- zscoreMatrix(humanBC, joincount.result) zscorePlot(matrix)