--- title: "FuseSOM package manual" author: "Elijah Willie" date: "`r Sys.Date()`" output: BiocStyle::html_document: toc: true vignette: > %\VignetteIndexEntry{FuseSOM package manual} %\VignetteEngine{knitr::rmarkdown} \usepackage[utf8]{inputenc} --- ```{r knitr-options, echo=FALSE, message=FALSE, warning=FALSE} library(knitr) opts_chunk$set(fig.align = 'center', fig.width = 6, fig.height = 5, dev = 'png') ``` # Installation ```{r, eval = FALSE} if (!require("BiocManager")) install.packages("BiocManager") BiocManager::install("FuseSOM") ``` # Introduction A correlation based multiview self organizing map for the characterization of cell types (`FuseSOM`) is a tool for unsupervised clustering. `FuseSOM` is robust and achieves high accuracy by combining a `Self Organizing Map` architecture and a `Multiview` integration of correlation based metrics to cluster highly multiplexed in situ imaging cytometry assays. The `FuseSOM` pipeline has been streamlined and accepts currently used data structures including `SingleCellExperiment` and `SpatialExperiment` objects as well as `DataFrames`. # Disclaimer This is purely a tool generated for clustering and as such it does not provide any means for QC and feature selection. It is advisable that the user first use other tools for quality control and feature selection before running `FuseSOM`. # Getting Started ## `FuseSOM` Matrix Input If you have a matrix containing expression data that was QCed and normalised by some other tool, the next step is to run the `FuseSOM` algorithm.This can be done by calling the `runFuseSOM()` function which takes in the matrix of interest where the columns are markers and the rows are observations, the makers of interest (if this is not provided, it is assumed that all columns are markers), and the number of clusters. ```{r, message=FALSE, warning=FALSE} # load FuseSOM library(FuseSOM) ``` Next we will load in the [`Risom et al`](https://www.sciencedirect.com/science/article/pii/S0092867421014860?via%3Dihub) dataset and run it through the FuseSOM pipeline. This dataset profiles the spatial landscape of ductal carcinoma in situ (DCIS), which is a pre-invasive lesion that is thought to be a precursor to invasive breast cancer (IBC). The key conclusion of this manuscript (amongst others) is that spatial information about cells can be used to predict disease progression in patients.We will also be using the markers used in the original study. ```{r} # load in the data data("risom_dat") # define the markers of interest risomMarkers <- c('CD45','SMA','CK7','CK5','VIM','CD31','PanKRT','ECAD', 'Tryptase','MPO','CD20','CD3','CD8','CD4','CD14','CD68','FAP', 'CD36','CD11c','HLADRDPDQ','P63','CD44') # we will be using the manual_gating_phenotype as the true cell type to gauge # performance names(risom_dat)[names(risom_dat) == 'manual_gating_phenotype'] <- 'CellType' ``` Now that we have loaded the data and define the markers of interest. We can run the `FuseSOM` algorithm. We have provided a function `runFuseSOM` that runs the pipeline from top to bottom and returns the cluster labels as well as the `Self Organizing Map` model. ```{r} risomRes <- runFuseSOM(data = risom_dat, markers = risomMarkers, numClusters = 23) ``` Lets look at the distribution of the clusters. ```{r} # get the distribution of the clusters table(risomRes$clusters)/sum(table(risomRes$clusters)) ``` Looks like `cluster_1` has about $32\%$ of the cells which is interesting. Next, lets generate a heatmap of the marker expression for each cluster. ```{r} risomHeat <- FuseSOM::markerHeatmap(data = risom_dat, markers = risomMarkers, clusters = risomRes$clusters, clusterMarkers = TRUE) ``` ## Using `FuseSOM` to estimate the number of clusters `FuseSOM` also provides functionality for estimating the number of clusters in a dataset using three classes of methods including: 1. Discriminant based method. + A method developed in house based on discriminant based maximum clusterability projection pursuit 2. Distance based methods which includes: + The Gap Statistic + The Jump Statistic + The Slope Statistic + The Within Cluster Dissimilarity Statistic + The Silhouette Statistic We can estimate the number of clusters using the `estimateNumCluster`. Run `help(estimateNumCluster)` to see it's complete functionality. ```{r} # lets estimate the number of clusters using all the methods # original clustering has 23 clusters so we will set kseq from 2:25 # we pass it the som model generated in the previous step risomKest <- estimateNumCluster(data = risomRes$model, kSeq = 2:25, method = c("Discriminant", "Distance")) ``` We can then use this result to determine the best number of clusters for this dataset based on the different metrics. The `FuseSOM` package provides a plotting function (`optiPlot`) which generates an elbow plot with the optimal value for the number of clusters for the distance based methods. See below ```{r} # what is the best number of clusters determined by the discriminant method? # optimal number of clusters according to the discriminant method is 7 risomKest$Discriminant # we can plot the results using the optiplot function pSlope <- optiPlot(risomKest, method = 'slope') pSlope pJump <- optiPlot(risomKest, method = 'jump') pJump pWcd <- optiPlot(risomKest, method = 'wcd') pWcd pGap <- optiPlot(risomKest, method = 'gap') pGap pSil <- optiPlot(risomKest, method = 'silhouette') pSil ``` From the plots, we see that the `Jump` statistics almost perfectly capture the number of clusters. The `Gap` method is a close second with $15$ clusters. All the other methods significantly underestimates the number of clusters. ## `FuseSOM` Sinlge Cell Epxeriment object as input. The `FuseSOM` algorithm is also equipped to take in a `SingleCellExperiment` object as input. The results of the pipeline will be written to either the metada or the colData fields. See below. First we create a `SingleCellExperiment` object ```{r, message=FALSE, warning=FALSE} library(SingleCellExperiment) # create a singelcellexperiment object colDat <- risom_dat[, setdiff(colnames(risom_dat), risomMarkers)] sce <- SingleCellExperiment(assays = list(counts = t(risom_dat)), colData = colDat) sce ``` Next we pass it to the `runFuseSOM()` function. Here, we can provide the assay in which the data is stored and what name to store the clusters under in the colData section. Note that the `Self Organizing Map` that is generated will be stored in the metadata field. ```{r} risomRessce <- runFuseSOM(sce, markers = risomMarkers, assay = 'counts', numClusters = 23, verbose = FALSE) colnames(colData(risomRessce)) names(metadata(risomRessce)) ``` Notice how the there is now a clusters column in the colData and SOM field in the metadata. You can run this function again with a new set of cluster number. If you provide a new name for the clusters, it will be stored under that new column, else, it will overwrite the current clusters column. Running it again on the same object will overwrite the SOM field in the metadata. Just like before, lets plot the heatmap of the resulting clusters across all markers. ```{r} data <- risom_dat[, risomMarkers] # get the original data used clusters <- colData(risomRessce)$clusters # extract the clusters from the sce object # generate the heatmap risomHeatsce <- markerHeatmap(data = risom_dat, markers = risomMarkers, clusters = clusters, clusterMarkers = TRUE) ``` ## Using `FuseSOM` to estimate the number of clusters for single cell experiment objects Just like before, we can estimate the number of clusters ```{r} # lets estimate the number of clusters using all the methods # original clustering has 23 clusters so we will set kseq from 2:25 # now we pass it a singlecellexperiment object instead of the som model as before # this will return a singelcellexperiment object where the metatdata contains the # cluster estimation information risomRessce <- estimateNumCluster(data = risomRessce, kSeq = 2:25, method = c("Discriminant", "Distance")) names(metadata(risomRessce)) ``` Notice how the metadata now contains a `clusterEstimation` field which holds the results from the `estimateNumCluster()` function We can assess the results in a similar fashion as before ```{r} # what is the best number of clusters determined by the discriminant method? # optimal number of clusters according to the discriminant method is 8 metadata(risomRessce)$clusterEstimation$Discriminant # we can plot the results using the optiplot function pSlope <- optiPlot(risomRessce, method = 'slope') pSlope pJump <- optiPlot(risomRessce, method = 'jump') pJump pWcd <- optiPlot(risomRessce, method = 'wcd') pWcd pGap <- optiPlot(risomRessce, method = 'gap') pGap pSil <- optiPlot(risomRessce, method = 'silhouette') pSil ``` Again, we see that the `Jump` statistics almost perfectly capture the number of clusters. The `Gap` method is a close second with $15$ clusters. All the other methods significantly underestimates the number of clusters. ## `FuseSOM` Spatial Epxeriment object as input. The methodology for `Spatial Epxeriment` is exactly the same as that of `Single Cell Epxeriment` # sessionInfo() ```{r} sessionInfo() ```