View source: R/genesignatures.R
sig_heatmap | R Documentation |
Plot a heatmap for the selected gene signature on the provided data, with the possibility to compactly display also DE only genes
sig_heatmap(
vst_data,
my_signature,
res_data = NULL,
FDR = 0.05,
de_only = FALSE,
annovec,
title = "",
cluster_rows = TRUE,
cluster_cols = FALSE,
anno_colData = NULL,
center_mean = TRUE,
scale_row = FALSE
)
vst_data |
A |
my_signature |
A character vector, usually named, containing the genes which compose the gene signature |
res_data |
A |
FDR |
Numeric value between 0 and 1, the False Discovery Rate |
de_only |
Logical, whether to display only DE genes belonging to the pathway - defaults to FALSE |
annovec |
A named character vector, with the corresponding annotation across IDs |
title |
Character, title for the heatmap |
cluster_rows |
Logical, whether to cluster rows - defaults to TRUE |
cluster_cols |
Logical, whether to cluster column - defaults to FALSE. Recommended to be set to TRUE if de_only is also set to TRUE |
anno_colData |
Character vector, specifying the elements of the colData information to be displayed as a decoration of the heatmap. Can be a vector of any length, as long as these names are included as colData. Defaults to NULL, which would plot no annotation on the samples. |
center_mean |
Logical, whether to perform mean centering on the expression values. Defaults to TRUE, as it improves the general readability of the heatmap |
scale_row |
Logical, whether to perform row-based standardization of the expression values |
A plot based on the pheatmap
function
# with the well known airway package...
library("airway")
data("airway", package = "airway")
airway
dds_airway <- DESeq2::DESeqDataSetFromMatrix(assay(airway),
colData = colData(airway),
design = ~ cell + dex
)
## Not run:
dds_airway <- DESeq2::DESeq(dds_airway)
res_airway <- DESeq2::results(dds_airway)
vst_airway <- DESeq2::vst(dds_airway)
library(org.Hs.eg.db)
annovec <- mapIds(org.Hs.eg.db, rownames(dds_airway), "ENTREZID", "ENSEMBL")
mysignatures <- read_gmt(
"http://data.wikipathways.org/20190210/gmt/wikipathways-20190210-gmt-Homo_sapiens.gmt"
)
mysignature_name <- "Lung fibrosis%WikiPathways_20190210%WP3624%Homo sapiens"
library(pheatmap)
sig_heatmap(vst_airway,
mysignatures[[mysignature_name]],
res_data = res_airway,
de_only = TRUE,
annovec = annovec,
title = mysignature_name,
cluster_cols = TRUE
)
## End(Not run)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.