R/markerDetection-methods.R
makeDetectionMatrices.Rd
Create matrices of detection events for individual gene features or combination thereof (i.e., signatures).
makeMarkerDetectionMatrix(se, markers, threshold = 0, assay.type = "counts") makeMarkerProportionScree(matrix) makeSignatureDetectionMatrix(matrix, object) makeMarkerProportionMatrix(se, cluster.col, assay.type = "counts", threshold = 0)
se | An object of class inheriting from |
---|---|
markers | A character vector, subset of |
threshold | Value above which the marker is considered detected. |
assay.type | A string specifying which assay values to use, e.g., |
matrix | A logical matrix indicating the presence of each marker (row) in each sample (column). |
object | A collection of signatures inheriting from " |
cluster.col | Name of a column in |
makeMarkerDetectionMatrix
A logical
matrix indicating detectable levels of each marker in each sample.
makeSignatureDetectionMatrix
A logical
matrix indicating detectable levels of each signature in each sample.
makeMarkerProportionMatrix
A matrix indicating for each feature the proportion of samples expressing detectable levels in each cluster.
makeMarkerProportionScree
A double
vector indicating the proportion of samples positive for markers in the input logical
matrix.
The makeMarkerDetectionMatrix
function returns a marker by sample matrix declaring a feature as detected if it is detected above a given threshold in a specific assay (e.g., "counts"
, "logcounts"
).
The makeSignatureDetectionMatrix
function returns a sample by signature matrix declaring a signature (composed of one or more gene features) as detected if all the associated features are detected.
The makeMarkerProportionMatrix
function returns a signature by cluster matrix indicating the proportion of samples expressing detectable levels of each marker in individual predefined clusters.
The makeMarkerProportionScree
function compute the 'cumulative' (i.e., combined) detection rate of markers:
the proportion of samples with detectable levels of the first marker, both of the first two markers, etc.
# Example data ---- library(SummarizedExperiment) nsamples <- 100 u <- matrix(rpois(20000, 2), ncol=nsamples) rownames(u) <- paste0("Gene", sprintf("%03d", seq_len(nrow(u)))) se <- SummarizedExperiment(assays=list(counts=u)) geneLists <- list( "Cell type 1" = c("Gene001", "Gene002"), "Cell type 2" = c("Gene003", "Gene004") ) bs <- as(geneLists, "Sets") # Example usage ---- markerMatrix <- makeMarkerDetectionMatrix(se, ids(elementInfo(bs))) signatureMatrix <- makeSignatureDetectionMatrix(markerMatrix, bs) tab <- makeMarkerProportionScree(markerMatrix) plot(tab$markers, tab$proportion, ylim=c(0, 1), main="Combined detection scree")se$cluster <- factor(sample(c("A", "B"), ncol(se), TRUE)) proportionMatrix <- makeMarkerProportionMatrix( se, cluster.col="cluster", assay.type="counts", threshold=0) head(proportionMatrix)#> cluster #> feature A B #> Gene001 0.8035714 0.9545455 #> Gene002 0.8571429 0.7954545 #> Gene003 0.9107143 0.8863636 #> Gene004 0.8928571 0.7727273 #> Gene005 0.7857143 0.9090909 #> Gene006 0.8928571 0.8863636