vsclust 1.0.0
Clustering is a method to identify common pattern in highly dimensional data. This can be for example genes or proteins with similar quantitative changes, thus providing insights into the affected biological pathways.
Despite of numerous clustering algorithms, they do not account for feature variance, i.e. the uncertainty in the measurements across the different experimental conditions. VSClust determines the characteristic patterns in high-dimensional data while accounting for feature variance that is given through replicated measurements.
Here, we present an example script to run the full clustering analysis using
the vsclust
library. The same can be done by running the Shiny app (e.g. via
its docker image or on ), or the
corresponding command line script. For the source code, see
.
Use the common Bioconductor commands for installation:
if (!require("BiocManager", quietly = TRUE))
install.packages("BiocManager")
BiocManager::install("vsclust")
The full functionality of this vignette can be obtained by additionally
installing and loading the packages matrixStats
and clusterProfiler
Here, we define the different parameters for the example data set
protein_expressions
. In the command-line version of VSClust (“runVSClust.R”),
they can be given via yaml file.
Comments:
A. Data sets with different numbers of replicates per condition need to be adapted to contain the same number of columns per condition. These can be done by either removing excess replicates or adding empty columns.
B. We assume the input data to be of the following format: A1, B1, C1, …, A2, B2, C2, …, where letters denote sample type and numbers are the different replicates.
C. If you prefer to estimate feature variance different, use averages and add
an estimate for the standard deviation as last column. You will need to set the
last option of PreparedForVSClust
to FALSE
.
D. If you don’t have replicates, use the same format as in C. and set the standard deviations to 1.
#### Input parameters, only read when now parameter file was provided
## All principal parameters for running VSClust can be defined as in the
## shinyapp at computproteomics.bmb.sdu.dk/Apps/VSClust
# name of study
Experiment <- "ProtExample"
# Number of replicates/sample per different experimental condition (sample
# type)
NumReps <- 3
# Number of different experimental conditions (e.g. time points or sample
# types)
NumCond <- 4
# Paired or unpaired statistical tests when carrying out LIMMA for
# statistical testing
isPaired <- FALSE
# Number of threads to accelerate the calculation (use 1 in doubt)
cores <- 1
# If 0 (default), then automatically estimate the cluster number for the
# vsclust
# run from the Minimum Centroid Distance
PreSetNumClustVSClust <- 0
# If 0 (default), then automatically estimate the cluster number for the
# original fuzzy c-means from the Minimum Centroid Distance
PreSetNumClustStand <- 0
# max. number of clusters when estimating the number of clusters. Higher
# numbers can drastically extend the computation time.
maxClust <- 10
At first, we load the example proteomics data set and carry out statistical
testing of all conditions version the first based on the LIMMA moderated t-test.
The data consists of mice fed with four different diets (high fat, TTA, fish oil
and TTA\(+\)fish oil).
Understand more about the data set with ?protein_expressions
This will calculate the false discovery rates for the differentially regulated features (pairwise comparisons versus the first “high fat” condition) and most importantly, their expected individual variances, to be used in the variance-sensitive clustering. These variances can also be uploaded separately via a last column containing them as individual standard deviations.
The PrepareForVSClust
function also creates a PCA plot to assess variability
and control whether the samples have been loaded correctly (replicated samples
should form groups).
After estimating the standard deviations, the matrix consists of the averaged quantitative feature values and a last column for the standard deviations of the features.
data(protein_expressions)
dat <- protein_expressions
#### running statistical analysis and estimation of individual variances
statOut <- PrepareForVSClust(dat, NumReps, NumCond, isPaired, TRUE)
dat <- statOut$dat
Sds <- dat[,ncol(dat)]
cat(paste("Features:",nrow(dat),"\nMissing values:",
sum(is.na(dat)),"\nMedian standard deviations:",
round(median(Sds,na.rm=TRUE),digits=3)))
## Features: 574
## Missing values: 0
## Median standard deviations: 0.22
## Write output into file
write.csv(statOut$statFileOut,
paste("",Experiment,"statFileOut.csv",sep=""))
There is no simple way to find the optimal number of clusters in a data set. For obtaining this number, we run the clustering for different cluster numbers and evaluate them via so-called validity indices, which provide information about suitable cluster numbers. VSClust uses mainly the “Maximum centroid distances” that denotes the shortest distance between any of the centroids. Alternatively, one can inspect the Xie Beni index.
The output of estimClustNum
contains the suggestion for the number of clusters.
We further visualize the outcome.
#### Estimate number of clusters with maxClust as maximum number clusters
#### to run the estimation with
ClustInd <- estimClustNum(dat, maxClust, cores)
## Running cluster number3
## Running cluster number4
## Running cluster number5
## Running cluster number6
## Running cluster number7
## Running cluster number8
## Running cluster number9
## Running cluster number10
#### Use estimate cluster number or use own
if (PreSetNumClustVSClust == 0)
PreSetNumClustVSClust <- optimalClustNum(ClustInd)
if (PreSetNumClustStand == 0)
PreSetNumClustStand <- optimalClustNum(ClustInd, method="FCM")
#### Visualize
estimClust.plot(ClustInd)
Now we run the clustering again with the optimal parameters from the estimation. One can take alternative numbers of clusters corresponding to large decays in the Minimum Centroid Distance or low values of the Xie Beni index.
First, we carry out the variance-sensitive method
#### Run clustering (VSClust and standard fcm clustering
ClustOut <- runClustWrapper(dat,
PreSetNumClustVSClust,
NULL,
VSClust=TRUE,
cores)
Bestcl <- ClustOut$Bestcl
VSClust_cl <- Bestcl
#ClustOut$p
## Write clustering results (VSClust)
write.csv(data.frame(cluster=Bestcl$cluster,
ClustOut$outFileClust,
isClusterMember=rowMaxs(Bestcl$membership)>0.5,
maxMembership=rowMaxs(Bestcl$membership),
Bestcl$membership),
paste(Experiment,
"FCMVarMResults",
Sys.Date(),
".csv",
sep=""))
## Write coordinates of cluster centroids
write.csv(Bestcl$centers,
paste(Experiment,
"FCMVarMResultsCentroids",
Sys.Date(),
".csv",
sep=""))
We see that most of the difference are between TTA diets and the rest. This shows that the TTA fatty acids have strong impact on the organisms. Cluster three shows the proteins that a commonly lower abundant in mice fed with fish oil and thus are related to biological processes affected this particular diet.
For comparison, this is the clustering using standard fuzzy c-means of the means over the replicates.
ClustOut <- runClustWrapper(dat, PreSetNumClustStand, NULL, VSClust=FALSE,
cores)
Bestcl <- ClustOut$Bestcl
## Write clustering results (standard fcm)
write.csv(data.frame(cluster=Bestcl$cluster,
ClustOut$outFileClust,
isClusterMember=rowMaxs(Bestcl$membership)>0.5,
maxMembership=rowMaxs(Bestcl$membership),
Bestcl$membership),
paste(Experiment,
"FCMResults",
Sys.Date(),
".csv",
sep=""))
## Write coordinates of cluster centroids
write.csv(Bestcl$centers, paste(Experiment,
"FCMResultsCentroids",
Sys.Date(),
".csv",
sep=""))
Here, the clusters look rather similar. VSClust best performs for larger numbers of different experimental conditions (one finds major improvements for \(D>6\)). For a 4-dimensional data set, the algorithm mostly filters out features with very high variance levels, making them unsuitable for belonging to a particular cluster.
This analysis is then followed by evaluating the features (here proteins) of
each cluster for their biological relevance. This can be done by functional
analysis with e.g. the clusterProfiler
package.
sessionInfo()
## R version 4.2.1 (2022-06-23)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 20.04.5 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.16-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.16-bioc/R/lib/libRlapack.so
##
## locale:
## [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
## [3] LC_TIME=en_GB LC_COLLATE=C
## [5] LC_MONETARY=en_US.UTF-8 LC_MESSAGES=en_US.UTF-8
## [7] LC_PAPER=en_US.UTF-8 LC_NAME=C
## [9] LC_ADDRESS=C LC_TELEPHONE=C
## [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C
##
## attached base packages:
## [1] stats4 stats graphics grDevices utils datasets methods
## [8] base
##
## other attached packages:
## [1] clusterProfiler_4.6.0 MultiAssayExperiment_1.24.0
## [3] SummarizedExperiment_1.28.0 Biobase_2.58.0
## [5] GenomicRanges_1.50.0 GenomeInfoDb_1.34.0
## [7] IRanges_2.32.0 S4Vectors_0.36.0
## [9] BiocGenerics_0.44.0 MatrixGenerics_1.10.0
## [11] matrixStats_0.62.0 vsclust_1.0.0
## [13] BiocStyle_2.26.0
##
## loaded via a namespace (and not attached):
## [1] ggtree_3.6.0 fgsea_1.24.0 colorspace_2.0-3
## [4] gson_0.0.9 ellipsis_0.3.2 qvalue_2.30.0
## [7] XVector_0.38.0 aplot_0.1.8 farver_2.1.1
## [10] graphlayouts_0.8.3 ggrepel_0.9.1 bit64_4.0.5
## [13] scatterpie_0.1.8 AnnotationDbi_1.60.0 fansi_1.0.3
## [16] codetools_0.2-18 splines_4.2.1 cachem_1.0.6
## [19] GOSemSim_2.24.0 knitr_1.40 polyclip_1.10-4
## [22] jsonlite_1.8.3 GO.db_3.16.0 png_0.1-7
## [25] ggforce_0.4.1 shiny_1.7.3 BiocManager_1.30.19
## [28] compiler_4.2.1 httr_1.4.4 lazyeval_0.2.2
## [31] assertthat_0.2.1 Matrix_1.5-1 fastmap_1.1.0
## [34] limma_3.54.0 cli_3.4.1 later_1.3.0
## [37] tweenr_2.0.2 htmltools_0.5.3 tools_4.2.1
## [40] igraph_1.3.5 gtable_0.3.1 glue_1.6.2
## [43] GenomeInfoDbData_1.2.9 reshape2_1.4.4 dplyr_1.0.10
## [46] fastmatch_1.1-3 Rcpp_1.0.9 enrichplot_1.18.0
## [49] jquerylib_0.1.4 vctrs_0.5.0 Biostrings_2.66.0
## [52] nlme_3.1-160 ape_5.6-2 ggraph_2.1.0
## [55] xfun_0.34 stringr_1.4.1 mime_0.12
## [58] lifecycle_1.0.3 DOSE_3.24.0 zlibbioc_1.44.0
## [61] MASS_7.3-58.1 scales_1.2.1 tidygraph_1.2.2
## [64] promises_1.2.0.1 parallel_4.2.1 RColorBrewer_1.1-3
## [67] yaml_2.3.6 gridExtra_2.3 memoise_2.0.1
## [70] ggfun_0.0.7 ggplot2_3.3.6 downloader_0.4
## [73] HDO.db_0.99.1 yulab.utils_0.0.5 sass_0.4.2
## [76] stringi_1.7.8 RSQLite_2.2.18 highr_0.9
## [79] tidytree_0.4.1 BiocParallel_1.32.0 rlang_1.0.6
## [82] pkgconfig_2.0.3 bitops_1.0-7 evaluate_0.17
## [85] lattice_0.20-45 purrr_0.3.5 treeio_1.22.0
## [88] patchwork_1.1.2 shadowtext_0.1.2 cowplot_1.1.1
## [91] bit_4.0.4 tidyselect_1.2.0 plyr_1.8.7
## [94] magrittr_2.0.3 bookdown_0.29 R6_2.5.1
## [97] magick_2.7.3 generics_0.1.3 DelayedArray_0.24.0
## [100] DBI_1.1.3 withr_2.5.0 pillar_1.8.1
## [103] KEGGREST_1.38.0 RCurl_1.98-1.9 tibble_3.1.8
## [106] crayon_1.5.2 utf8_1.2.2 rmarkdown_2.17
## [109] viridis_0.6.2 grid_4.2.1 data.table_1.14.4
## [112] blob_1.2.3 digest_0.6.30 xtable_1.8-4
## [115] tidyr_1.2.1 httpuv_1.6.6 gridGraphics_0.5-1
## [118] BiocBaseUtils_1.0.0 munsell_0.5.0 ggplotify_0.1.0
## [121] viridisLite_0.4.1 bslib_0.4.0