Skip to contents

Get ICP cell cluster probability table(s)

Usage

GetCellClusterProbability.SingleCellExperiment(
  object,
  icp.run,
  icp.round,
  concatenate
)

# S4 method for class 'SingleCellExperiment'
GetCellClusterProbability(
  object,
  icp.run = NULL,
  icp.round = NULL,
  concatenate = TRUE
)

Arguments

object

An object of SingleCellExperiment class with ICP cell cluster probability tables saved in metadata(object)$coralysis$joint.probability. After running RunParallelDivisiveICP.

icp.run

ICP run(s) to retrieve from metadata(object)$coralysis$joint.probability. By default NULL, i.e., all are retrieved. Specify a numeric vector to retrieve a specific set of tables.

icp.round

ICP round(s) to retrieve from metadata(object)$coralysis$joint.probability. By default NULL, i.e., all are retrieved.

concatenate

Concatenate list of ICP cell cluster probability tables retrieved. By default TRUE, i.e., the list of ICP cell cluster probability tables is concatenated.

Value

A list with ICP cell cluster probability tables or a matrix with concatenated tables.

Examples

# Import package
suppressPackageStartupMessages(library("SingleCellExperiment"))

# Create toy SCE data
batches <- c("b1", "b2")
set.seed(239)
batch <- sample(x = batches, size = nrow(iris), replace = TRUE)
sce <- SingleCellExperiment(assays = list(logcounts = t(iris[,1:4])),  
                            colData = DataFrame("Species" = iris$Species, 
                                               "Batch" = batch))
colnames(sce) <- paste0("samp", 1:ncol(sce))

# Prepare SCE object for analysis
sce <- PrepareData(sce)
#> Converting object of `matrix` class into `dgCMatrix`. Please note that Coralysis has been designed to work with sparse data, i.e. data with a high proportion of zero values! Dense data will likely increase run time and memory usage drastically!
#> 4/4 features remain after filtering features with only zero values.

# Multi-level integration (just for highlighting purposes; use default parameters)
set.seed(123)
sce <- RunParallelDivisiveICP(object = sce, batch.label = "Batch", 
                              k = 2, L = 25, C = 1, train.k.nn = 10, 
                              train.k.nn.prop = NULL, use.cluster.seed = FALSE,
                              build.train.set = FALSE, ari.cutoff = 0.1, 
                             threads = 2)
#> 
#> Initializing divisive ICP clustering...
#> 
  |                                                                            
  |                                                                      |   0%
  |                                                                            
  |===                                                                   |   4%
  |                                                                            
  |======                                                                |   8%
  |                                                                            
  |=========                                                             |  12%
  |                                                                            
  |============                                                          |  17%
  |                                                                            
  |===============                                                       |  21%
  |                                                                            
  |==================                                                    |  25%
  |                                                                            
  |====================                                                  |  29%
  |                                                                            
  |=======================                                               |  33%
  |                                                                            
  |==========================                                            |  38%
  |                                                                            
  |=============================                                         |  42%
  |                                                                            
  |================================                                      |  46%
  |                                                                            
  |===================================                                   |  50%
  |                                                                            
  |======================================                                |  54%
  |                                                                            
  |=========================================                             |  58%
  |                                                                            
  |============================================                          |  62%
  |                                                                            
  |===============================================                       |  67%
  |                                                                            
  |==================================================                    |  71%
  |                                                                            
  |====================================================                  |  75%
  |                                                                            
  |=======================================================               |  79%
  |                                                                            
  |==========================================================            |  83%
  |                                                                            
  |=============================================================         |  88%
  |                                                                            
  |================================================================      |  92%
  |                                                                            
  |===================================================================   |  96%
  |                                                                            
  |======================================================================| 100%
#> 
#> Divisive ICP clustering completed successfully.
#> 
#> Predicting cell cluster probabilities using ICP models...
#> Prediction of cell cluster probabilities completed successfully.
#> 
#> Multi-level integration completed successfully.

# Get cluster probability for all ICP runs 
probs <- GetCellClusterProbability(object = sce, icp.round = 1, concatenate = TRUE) 
probs[1:10, 1:5]
#>                1          2         1          2          1
#> samp1  0.9643473 0.03565274 0.9638023 0.03619770 0.03965366
#> samp2  0.9264623 0.07353774 0.9255917 0.07440832 0.08002244
#> samp3  0.9519228 0.04807721 0.9512637 0.04873629 0.05290241
#> samp4  0.9268923 0.07310768 0.9260116 0.07398839 0.07971891
#> samp5  0.9692441 0.03075587 0.9687533 0.03124672 0.03435248
#> samp6  0.9698052 0.03019485 0.9693037 0.03069632 0.03394382
#> samp7  0.9587039 0.04129606 0.9581008 0.04189923 0.04573415
#> samp8  0.9524988 0.04750117 0.9518274 0.04817258 0.05249006
#> samp9  0.9153487 0.08465131 0.9144061 0.08559393 0.09170199
#> samp10 0.9268923 0.07310768 0.9260116 0.07398839 0.07971891