Title: | Single Cell Transcriptomics-Level Cytokine Activity Prediction and Estimation |
---|---|
Description: | Generates cell-level cytokine activity estimates using relevant information from gene sets constructed with the 'CytoSig' and the 'Reactome' databases and scored using the modified 'Variance-adjusted Mahalanobis (VAM)' framework for single-cell RNA-sequencing (scRNA-seq) data. 'CytoSig' database is described in: Jiang at al., (2021) <doi:10.1038/s41592-021-01274-5>. 'Reactome' database is described in: Gillespie et al., (2021) <doi:10.1093/nar/gkab1028>. The 'VAM' method is outlined in: Frost (2020) <doi:10.1093/nar/gkaa582>. |
Authors: | H. Robert Frost [aut], Azka Javaid [aut, cre] |
Maintainer: | Azka Javaid <[email protected]> |
License: | GPL (>= 2) |
Version: | 0.1.0 |
Built: | 2024-12-12 07:14:29 UTC |
Source: | CRAN |
Returns the gene set associated with input cytokine(s) from the CytoSig database given manually specified cytokine-specific output csv file(s) under the extdata directory with file name beginning with the specified cytokine (i.e., 'IL6_output.csv') as currently provided for the IL6 cytokine.
genesetCytoSig(cytokine.eval, file.name)
genesetCytoSig(cytokine.eval, file.name)
cytokine.eval |
Cytokine(s) associated with the query. |
file.name |
List of XML file(s) associated with the cytokine(s) beginning with the specific cytokine name. |
Dataframe consisting of genes and the associated log2 fold change values associated with the specific cytokine(s).
file.name.cytosig1 <- system.file("extdata", "IL6_output.csv", package = "scaper") genesetCytoSig(cytokine.eval = "IL6", file.name = file.name.cytosig1) %>% head(10)
file.name.cytosig1 <- system.file("extdata", "IL6_output.csv", package = "scaper") genesetCytoSig(cytokine.eval = "IL6", file.name = file.name.cytosig1) %>% head(10)
Returns the gene set consisting of genes on the Reactome pathway hierarchy given input cytokine(s) and specified list of pathway hierarchy xml file(s) with file name beginning with the specified cytokine (i.e., 'IL6_Interleukin6_signaling.xml') as currently provided for the IL6 cytokine.
genesetReactome(cytokine.eval, file.name)
genesetReactome(cytokine.eval, file.name)
cytokine.eval |
Cytokine(s) associated with the query. |
file.name |
List of XML file(s) associated with the cytokine(s) beginning with the specific cytokine name. |
Dataframe consisting of list of genes on the molecular pathway associated with the specified cytokine(s).
file.name.reactome1 <- system.file("extdata", "IL6_Interleukin6_signaling.xml", package = "scaper") genesetReactome(cytokine.eval = "IL6", file.name = file.name.reactome1) %>% head(10)
file.name.reactome1 <- system.file("extdata", "IL6_Interleukin6_signaling.xml", package = "scaper") genesetReactome(cytokine.eval = "IL6", file.name = file.name.reactome1) %>% head(10)
Computes cell-level estimates of cytokine activity for a normalized scRNA-seq
count matrix using the SCAPE method. SCAPE activity estimates are computed by scoring weighted genes
sets from the CytoSig or Reactome databases using the VAM::vamForCollection()
function. Individual gene sets for subsequent scoring can be reconstructed using
the genesetCytoSig
and the genesetReactome
functions
for the CytoSig and the Reactome database, respectively.
scape(counts.matrix, database = "cytosig", cytokine = "all")
scape(counts.matrix, database = "cytosig", cytokine = "all")
counts.matrix |
A |
database |
Database used for gene set construction and set scoring.
|
cytokine |
Vector of cytokine names to score for activity. The default value
of "all" will score all 41 cytokines supported by CytoSig or 31 supported by Reactome. Please see
function |
A matrix consisting of the cell-level cytokine activity scores for
cytokines.
genesetCytoSig
, genesetReactome
library(Seurat) library(SeuratObject) pbmc_small <- NormalizeData(pbmc_small) counts.matrix <- as.data.frame(t(as.matrix(pbmc_small@assays$RNA@data))) CytoSig.score.output <- scape(counts.matrix = counts.matrix, database = "cytosig") head(CytoSig.score.output)[,1:3] CytoSig.score.output.specific <- scape(counts.matrix = counts.matrix, database = "cytosig", cytokine = c("IL4", "IL13")) head(CytoSig.score.output.specific)
library(Seurat) library(SeuratObject) pbmc_small <- NormalizeData(pbmc_small) counts.matrix <- as.data.frame(t(as.matrix(pbmc_small@assays$RNA@data))) CytoSig.score.output <- scape(counts.matrix = counts.matrix, database = "cytosig") head(CytoSig.score.output)[,1:3] CytoSig.score.output.specific <- scape(counts.matrix = counts.matrix, database = "cytosig", cytokine = c("IL4", "IL13")) head(CytoSig.score.output.specific)
Computes cell-level estimates of cytokine activity for a scRNA-seq Seurat
count matrix using the scapeForSeurat method. SCAPE activity estimates are computed by scoring weighted gene
sets from the CytoSig or Reactome databases using the Variance-adjusted Mahalanobis (VAM) method
as implemented in the VAM::vamForSeurat()
function. Individual gene sets for subsequent scoring can be reconstructed using
the genesetCytoSig
and the genesetReactome
functions
for the CytoSig and the Reactome database, respectively.
scapeForSeurat( seurat.object, database = "cytosig", cytokine = "all", normalize = TRUE )
scapeForSeurat( seurat.object, database = "cytosig", cytokine = "all", normalize = TRUE )
seurat.object |
Seurat counts matrix. |
database |
Database used for gene set construction and set scoring.
|
cytokine |
Vector of cytokine names to score for activity. The default value of "all"
will score all 41 cytokines supported by CytoSig or 31 supported by Reactome. Please see
function |
normalize |
Boolean indicator for whether normalization should be performed before performing gene set scoring. |
Seurat object consisting of cell-level cytokine activity scores returned as a separate assay (scape for scoring via the CytoSig database and VAMcdf for scoring via the Reactome database).
genesetCytoSig
, genesetReactome
, scape
library(SeuratObject) CytoSig.score.output.all <- scapeForSeurat(seurat.object = pbmc_small, database = "cytosig", cytokine = "all", normalize=TRUE) (as.data.frame(CytoSig.score.output.all@assays$scape@data))[1:6,1:3] CytoSig.score.output.specific <- scapeForSeurat(seurat.object = pbmc_small, database = "cytosig", cytokine = c("IL4", "IL13"), normalize=TRUE) (as.data.frame(CytoSig.score.output.specific@assays$scape@data))[,1:3]
library(SeuratObject) CytoSig.score.output.all <- scapeForSeurat(seurat.object = pbmc_small, database = "cytosig", cytokine = "all", normalize=TRUE) (as.data.frame(CytoSig.score.output.all@assays$scape@data))[1:6,1:3] CytoSig.score.output.specific <- scapeForSeurat(seurat.object = pbmc_small, database = "cytosig", cytokine = c("IL4", "IL13"), normalize=TRUE) (as.data.frame(CytoSig.score.output.specific@assays$scape@data))[,1:3]
Returns the names of the cytokines supported by either the CytoSig or the Reactome databases.
supportedCytokines(database = "cytosig")
supportedCytokines(database = "cytosig")
database |
Database used for gene set construction and set scoring.
|
List of cytokines associated with the CytoSig or the Reactome databases.
supportedCytokines(database = "cytosig") supportedCytokines(database = "reactome")
supportedCytokines(database = "cytosig") supportedCytokines(database = "reactome")