specificityScore | R Documentation |
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 signature 'matrix'
specificityScore(
x,
method = c("tau", "TSI", "counts"),
group = NULL,
thresh = 0,
expr_values = "logcounts",
na.rm = FALSE
)
## S4 method for signature 'SummarizedExperiment'
specificityScore(
x,
method = c("tau", "TSI", "counts"),
group = NULL,
thresh = 0,
expr_values = "logcounts",
na.rm = FALSE
)
x |
Expression data, either a |
method |
|
group |
|
thresh |
|
expr_values |
Integer scalar or string indicating which assay of
|
na.rm |
Logical scalar. If |
A numeric
vector of length nrow(x)
with gene scores.
Michael Stadler
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
x <- rbind(g1 = runif(5),
g2 = c(1, 0, 0, 0, 0),
g3 = c(.6, .1, .1, .1, .1))
specificityScore(x)
specificityScore(x, method = "TSI")
specificityScore(x, method = "counts", thresh = 0.5)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.