Calculate expression specificity scores for genes that quantify specific expression of a gene in groups of samples (e.g. from different tissues).

specificityScore(
  x,
  method = c("tau", "TSI", "counts"),
  group = NULL,
  thresh = 0,
  expr_values = "logcounts",
  na.rm = FALSE
)

# S4 method for class 'matrix'
specificityScore(
  x,
  method = c("tau", "TSI", "counts"),
  group = NULL,
  thresh = 0,
  expr_values = "logcounts",
  na.rm = FALSE
)

# S4 method for class 'SummarizedExperiment'
specificityScore(
  x,
  method = c("tau", "TSI", "counts"),
  group = NULL,
  thresh = 0,
  expr_values = "logcounts",
  na.rm = FALSE
)

Arguments

x

Expression data, either a matrix with expression values for genes (rows) in each sample (columns), or a SummarizedExperiment or SingleCellExperiment object containing such expression data in one of the assays (selected by expr_values).

method

character scalar selecting the type of expression specificity score to be calculated. One of: "tau", "TSI", "counts". See "Details" for method-specific information.

group

character or factor of length ncol(x) that groups the measurements into clusters or tissues, for which expression specificity scores are to be calculated. If NULL (the default), each column of x is assumed to be its own group. If multiple columns belong to the same group, these columns are first aggregated (averaged) before score calculations.

thresh

numeric scalar defining the expression threshold. Values greater than this threshold are interpreted as expressed (used only for some of the methods, see "Details").

expr_values

Integer scalar or string indicating which assay of x contains the expression values, for x of type SummarizedExperiment or SingleCellExperiment. Ignored if x is a matrix.

na.rm

Logical scalar. If TRUE, NA values are excluded in the calculations.

Value

A numeric vector of length nrow(x) with gene scores.

References

For a review of tissue-specificity scores, see: "A benchmark of gene expression tissue-specificity metrics" Nadezda Kryuchkova-Mostacci and Marc Robinson-Rechavi Brief Bioinform. 2017 Mar; 18(2): 205–214. doi: 10.1093/bib/bbw008, PMCID: PMC5444245, PMID: 26891983

Author

Michael Stadler

Examples

x <- rbind(g1 = runif(5),
           g2 = c(1, 0, 0, 0, 0),
           g3 = c(.6, .1, .1, .1, .1))
specificityScore(x)
#> calculating 'tau' for 3 genes and 5 groups
#>        g1        g2        g3 
#> 0.4584552 1.0000000 0.8333333 
specificityScore(x, method = "TSI")
#> calculating 'TSI' for 3 genes and 5 groups
#>        g1        g2        g3 
#> 0.3158381 1.0000000 0.6000000 
specificityScore(x, method = "counts", thresh = 0.5)
#> calculating 'counts' for 3 genes and 5 groups
#> g1 g2 g3 
#>  3  1  1