To install this package, start R (version > “4.0”) and enter:
If you have any FEAST-related questions, please post to the GitHub Issue section of FEAST at https://github.com/suke18/FEAST/issues, which will be helpful for the construction of FEAST.
Cell clustering is one of the most important and commonly performed tasks in single-cell RNA sequencing (scRNA-seq) data analysis. An important step in cell clustering is to select a subset of genes (referred to as “features”), whose expression patterns will then be used for downstream clustering. A good set of features should include the ones that distinguish different cell types, and the quality of such set could have significant impact on the clustering accuracy. All existing scRNA-seq clustering tools include a feature selection step relying on some simple unsupervised feature selection methods, mostly based on the statistical moments of gene-wise expression distributions. In this work, we develop a feature selection algorithm named FEAture SelecTion (FEAST), which provides more representative features.
We use Yan
dataset(Yan et al.
2013) for the demonstration. This Yan
dataset
includes 6 cell types about human preimplantation embryos and embryonic
stem cells. The users can run the FEAST or
FEAST_fast function to obtain the gene rankings (from
the most significant to the least significant genes) of These featured
genes are deemed to better contribute to the clustering accuracy. Here,
we demonstrate that some of the top 9 selected genes are very
informative (or known as marker genes e.g., CYYR1 ).
library(FEAST)
#> Warning: multiple methods tables found for 'setequal'
#> Warning: replacing previous import 'BiocGenerics::setequal' by
#> 'S4Vectors::setequal' when loading 'IRanges'
#> Warning: replacing previous import 'BiocGenerics::setequal' by
#> 'S4Vectors::setequal' when loading 'GenomeInfoDb'
#> Warning: replacing previous import 'BiocGenerics::setequal' by
#> 'S4Vectors::setequal' when loading 'GenomicRanges'
#> Warning: replacing previous import 'BiocGenerics::setequal' by
#> 'S4Vectors::setequal' when loading 'XVector'
#> Warning: replacing previous import 'BiocGenerics::setequal' by
#> 'S4Vectors::setequal' when loading 'SummarizedExperiment'
#> Warning: replacing previous import 'BiocGenerics::setequal' by
#> 'S4Vectors::setequal' when loading 'S4Arrays'
#> Warning: replacing previous import 'BiocGenerics::setequal' by
#> 'S4Vectors::setequal' when loading 'DelayedArray'
#> Warning: replacing previous import 'BiocGenerics::setequal' by
#> 'S4Vectors::setequal' when loading 'SparseArray'
#> Warning: replacing previous import 'S4Arrays::read_block' by
#> 'DelayedArray::read_block' when loading 'SummarizedExperiment'
data(Yan)
k = length(unique(trueclass))
Y = process_Y(Y, thre = 2)
# The core function.
ixs = FEAST(Y, k=k)
# look at the features
Ynorm = Norm_Y(Y)
par(mfrow = c(3,3))
for (i in 1:9){
tmp_ix = ixs[i]
tmp_gene = rownames(Ynorm)[tmp_ix]
boxplot(as.numeric(Ynorm[tmp_ix, ])~trueclass, main = tmp_gene, xlab="", ylab="", las=2)
}
FEAST requires a count expression matrix with rows representing the
genes and columns representing the cells. To preprocess the count
matrix, FEAST filters out the genes depending on the dropout rates
(sometimes known as zero-expression rates). Alternatively, the users can
input the normalized matrix (usually log transformed of relative
abundance expression matrix (scaled by average sequencing depth)). Here,
we use one template Yan
(Yan et al.
2013) dataset. It is important to note that FEAST requires the
number of clusters (k) as an input, which can be
obtained by the prior knowledge of the cell types.
The Yan
dataset is loaded with two objects including the
count expression matrix (Y) and corresponding cell type information
(trueclass). Other sample datasets can be downloaded at https://drive.google.com/drive/u/0/folders/1SRT7mrX7ziJoSjuFLLkK8kjnUsJrabVM.
If the users want to use one Deng dataset(Deng et
al. 2014), the users can load the data and access the phenotype
information by using colData(Deng)
function, and the count
expression matrix by using assay(Deng,"counts")
function
from the SummarizedExperiment Bioconductor package.
To preprocess the count expression matrix (Y
), FEAST
filters out the genes based on the dropout rate (sometimes known as
zero-expression rate). This consensus clustering step will output the
initial clustering labels based on a modified CSPA algorithm. It will
find the cells that tightly clustered together and possibly filter out
the cells that are less correlated to the initial clustering
centers.
After the consensus clustering step, FEAST will generate the initial clustering labels with confidence. Based on the initial clustering outcomes, FEAST further infers the gene-level significance by using F-statistics followed by the ranking of the full list of genes.
After ranking the genes by F-statistics, FEAST curates a series of
top number of features (genes). Then, FEAST input these top (by default
top = 500, 1000, 2000
) features into some established
clustering algorithm such as SC3 (Kiselev et al.
2017), which is regarded as the most accurate scRNA-seq
clustering algorithm (Duò, Robinson, and Soneson
2018). It is worth noting that one needs to confirm the
k to be the same for every iteration.
After the clustering steps, with the case of unknown cell types, FEAST evaluates the quality of the feature set by computing the average distance between each cell and the clustering centers, which is the same as the MSE. It is noting that FEAST uses all the features (genes) for distance calculation for the purpose of fair comparisons. The clustering centers are obtained by using previous clustering step.
It is noted that FEAST is not limited to SC3 clustering. It also shows superior performance on other clustering methods such as TSCAN. The code for running SC3 and TSCAN are embedded in FEAST package, which can be access by SC3_Clust and TSCAN_Clust.
## clustering step
tops = c(500, 1000, 2000)
cluster_res = NULL
for (top in tops){
tmp_ixs = ixs[1:top]
tmp_markers = rownames(Y)[tmp_ixs]
tmp_res = TSCAN_Clust(Y, k = k, input_markers = tmp_markers)
#tmp_res = SC3_Clust(Y, k = k, input_markers = tmp_markers)
cluster_res[[toString(top)]] = tmp_res
}
## validation step
Ynorm = Norm_Y(Y)
mse_res = NULL
for (top in names(cluster_res)){
tmp_res = cluster_res[[top]]
tmp_cluster = tmp_res$cluster
tmp_mse = cal_MSE(Ynorm = Ynorm, cluster = tmp_cluster)
mse_res = c(mse_res, tmp_mse)
}
names(mse_res) = names(cluster_res)
After the validation step, the feature set associated with the smallest MSE value will be recommended for the further analysis. Here, we demonstrate that an improved clustering accuracy by specifying the optimal feature set. The benchmark comparison is the original setting of established clustering algorithm (e.g. SC3).
As demonstrated in the figure below, the first panel includes the PCA
illustration of the cell types. The second panel shows the clustering
result by the original SC3, in which some cells (inside the gray circle)
are mixed. The third panel shows the improved clustering outcomes by
inputing the optimized feature set into SC3. This is an example for the
Deng
dataset.
The users can run the Consensus
function and then run
the Select_Model_short_SC3
function; however,
Consensus
step could take a while especially for the large
dataset (sample size is greater than 2000).
data(Yan)
Y = process_Y(Y, thre = 2) # preprocess the data if needed
con_res = Consensus(Y, k=k)
mod_res = Select_Model_short_TSCAN(Y, cluster = con_res$cluster, top = c(200, 500, 1000, 2000))
# mod_res = Select_Model_short_SC3(Y, cluster = con_res$cluster, top = c(200, 500, 1000, 2000))
# to visualize the result, one needs to load ggpubr library.
library(ggpubr)
PP = Visual_Rslt(model_cv_res = mod_res, trueclass = trueclass)
print(PP$ggobj) # show the visualization plot.
The result figure shows that the MSE validation process is concordant with the clustering accuracy trend. It demonstrated that using top 1000 featured genes can result in the highest clustering accuracy with improvement than the original setting.
To quickly obtain the ranking orders of the features, the users can
directly apply FEAST
function, which works perfectly fun on
small dataset. For the large dataset, the users can consider change the
dimention reduction method to irlba by specifying
dim_reduce="irlba"
. Moreover, the users can resort to
FEAST_fast
function for fast calculation on the large
datasets. For extreme large dataset (sample size >5000), it will
split the dataset equally into chucks and apply FEAST algorithm in
parallel on individual splitted datasets. In this case, the users need
to specify two parameters split = TRUE
, and
batch_size = 1000
corresponding to the splitted size.
sessionInfo()
#> R version 4.4.2 (2024-10-31)
#> Platform: x86_64-pc-linux-gnu
#> Running under: Ubuntu 24.04.1 LTS
#>
#> Matrix products: default
#> BLAS: /usr/lib/x86_64-linux-gnu/openblas-pthread/libblas.so.3
#> LAPACK: /usr/lib/x86_64-linux-gnu/openblas-pthread/libopenblasp-r0.3.26.so; LAPACK version 3.12.0
#>
#> locale:
#> [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
#> [3] LC_TIME=en_US.UTF-8 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
#>
#> time zone: Etc/UTC
#> tzcode source: system (glibc)
#>
#> attached base packages:
#> [1] stats4 stats graphics grDevices utils datasets methods
#> [8] base
#>
#> other attached packages:
#> [1] ggpubr_0.6.0 ggplot2_3.5.1
#> [3] FEAST_1.15.0 SummarizedExperiment_1.37.0
#> [5] Biobase_2.67.0 GenomicRanges_1.59.0
#> [7] GenomeInfoDb_1.43.0 IRanges_2.41.0
#> [9] S4Vectors_0.45.0 BiocGenerics_0.53.1
#> [11] generics_0.1.3 MatrixGenerics_1.19.0
#> [13] matrixStats_1.4.1 BiocParallel_1.41.0
#> [15] mclust_6.1.1 BiocStyle_2.35.0
#>
#> loaded via a namespace (and not attached):
#> [1] bitops_1.0-9 rlang_1.1.4
#> [3] magrittr_2.0.3 e1071_1.7-16
#> [5] compiler_4.4.2 mgcv_1.9-1
#> [7] vctrs_0.6.5 combinat_0.0-8
#> [9] pkgconfig_2.0.3 crayon_1.5.3
#> [11] fastmap_1.2.0 backports_1.5.0
#> [13] XVector_0.47.0 labeling_0.4.3
#> [15] caTools_1.18.3 utf8_1.2.4
#> [17] promises_1.3.0 rmarkdown_2.29
#> [19] UCSC.utils_1.3.0 TSCAN_1.45.0
#> [21] purrr_1.0.2 xfun_0.49
#> [23] WriteXLS_6.7.0 zlibbioc_1.52.0
#> [25] cachem_1.1.0 jsonlite_1.8.9
#> [27] highr_0.11 later_1.3.2
#> [29] DelayedArray_0.33.1 broom_1.0.7
#> [31] irlba_2.3.5.1 parallel_4.4.2
#> [33] cluster_2.1.6 R6_2.5.1
#> [35] bslib_0.8.0 RColorBrewer_1.1-3
#> [37] car_3.1-3 rrcov_1.7-6
#> [39] jquerylib_0.1.4 Rcpp_1.0.13-1
#> [41] iterators_1.0.14 knitr_1.48
#> [43] splines_4.4.2 igraph_2.1.1
#> [45] httpuv_1.6.15 Matrix_1.7-1
#> [47] tidyselect_1.2.1 abind_1.4-8
#> [49] yaml_2.3.10 doParallel_1.0.17
#> [51] gplots_3.2.0 codetools_0.2-20
#> [53] doRNG_1.8.6 plyr_1.8.9
#> [55] lattice_0.22-6 tibble_3.2.1
#> [57] withr_3.0.2 shiny_1.9.1
#> [59] ROCR_1.0-11 evaluate_1.0.1
#> [61] proxy_0.4-27 TrajectoryUtils_1.15.0
#> [63] pillar_1.9.0 BiocManager_1.30.25
#> [65] carData_3.0-5 rngtools_1.5.2
#> [67] SC3_1.35.0 KernSmooth_2.23-24
#> [69] foreach_1.5.2 pcaPP_2.0-5
#> [71] munsell_0.5.1 scales_1.3.0
#> [73] fastICA_1.2-5.1 gtools_3.9.5
#> [75] xtable_1.8-4 class_7.3-22
#> [77] glue_1.8.0 pheatmap_1.0.12
#> [79] maketools_1.3.1 tools_4.4.2
#> [81] robustbase_0.99-4-1 sys_3.4.3
#> [83] ggsignif_0.6.4 buildtools_1.0.0
#> [85] mvtnorm_1.3-2 cowplot_1.1.3
#> [87] grid_4.4.2 tidyr_1.3.1
#> [89] colorspace_2.1-1 SingleCellExperiment_1.29.0
#> [91] nlme_3.1-166 GenomeInfoDbData_1.2.13
#> [93] Formula_1.2-5 cli_3.6.3
#> [95] fansi_1.0.6 S4Arrays_1.7.1
#> [97] dplyr_1.1.4 gtable_0.3.6
#> [99] DEoptimR_1.1-3 ggsci_3.2.0
#> [101] rstatix_0.7.2 sass_0.4.9
#> [103] digest_0.6.37 SparseArray_1.7.1
#> [105] farver_2.1.2 htmltools_0.5.8.1
#> [107] lifecycle_1.0.4 httr_1.4.7
#> [109] mime_0.12