Title: | Unsupervised Identification of Genes with Expression Loss in Subsets of Tumors |
---|---|
Description: | receptLoss identifies genes whose expression is lost in subsets of tumors relative to normal tissue. It is particularly well-suited in cases where the number of normal tissue samples is small, as the distribution of gene expression in normal tissue samples is approximated by a Gaussian. Originally designed for identifying nuclear hormone receptor expression loss but can be applied transcriptome wide as well. |
Authors: | Daniel Pique, John Greally, Jessica Mar |
Maintainer: | Daniel Pique <[email protected]> |
License: | GPL-3 + file LICENSE |
Version: | 1.19.0 |
Built: | 2024-11-27 05:10:58 UTC |
Source: | https://github.com/bioc/receptLoss |
This object contains a table of all known NHRs and was adapted from the 'guidetopharmacology' website (see references). It was joined with a bioMart table to include ensemble gene ids, which are commonly used gene symbols.
nhrs
nhrs
A tibble with 54 rows and 6 variables:
the HUGO gene nomenclature committee (HGNC) symbol (letters and numbers, ex. THRB)
the HUGO gene nomenclature committee (HGNC) symbol (a number, ex. 11799)
the HUGO gene nomenclature committee (HGNC) gene name (ex. "Thyroid hormone receptor beta")
the entrez gene id (a number, ex. 7068)
the ensembl gene id (ex. ENSG00000151090, always starts with ENSG)
words or gene symbols in the literature that refer to the same gene
http://www.guidetopharmacology.org/DATA/targets_and_families.csv
This function allows you to plot histograms of tumor and adj normal data
plotReceptLoss(exprMatrNml, exprMatrTum, rldf, geneName, addToTitle = "", clrs)
plotReceptLoss(exprMatrNml, exprMatrTum, rldf, geneName, addToTitle = "", clrs)
exprMatrNml |
A matrix of expression values from normal tissue. Each row is a gene, and each column is a patient or sample. Genes should be in same order as exprMatrTum. |
exprMatrTum |
A matrix of expression values from tumor tissue. Each row is a gene, and each column is a patient or sample. Genes should be in same order as exprMatrNml. |
rldf |
The dataframe output from running the receptLoss function |
geneName |
The name of the gene to plot. The name of the gene should correspond to a row name in both exprMatrNml and exprMatrTum matrices. |
addToTitle |
A string that can be added to the title, which includes the gene name. |
clrs |
Vector of length 2 containing colors to use for plot |
returns an object of class 'ggplot'
exprMatrNml <- matrix(abs(rnorm(100, mean=2)), nrow=10) exprMatrTum <- matrix(abs(rnorm(100)), nrow=10) geneNames <- paste0(letters[seq_len(nrow(exprMatrNml))], seq_len(nrow(exprMatrNml))) rownames(exprMatrNml) <- rownames(exprMatrTum) <- geneNames nSdBelow <- 2 minPropPerGroup <- .2 rl <- receptLoss(exprMatrNml, exprMatrTum, nSdBelow, minPropPerGroup) clrs <- c("#E78AC3", "#8DA0CB") plotReceptLoss(exprMatrNml, exprMatrTum, rl, geneName="g7", clrs=clrs)
exprMatrNml <- matrix(abs(rnorm(100, mean=2)), nrow=10) exprMatrTum <- matrix(abs(rnorm(100)), nrow=10) geneNames <- paste0(letters[seq_len(nrow(exprMatrNml))], seq_len(nrow(exprMatrNml))) rownames(exprMatrNml) <- rownames(exprMatrTum) <- geneNames nSdBelow <- 2 minPropPerGroup <- .2 rl <- receptLoss(exprMatrNml, exprMatrTum, nSdBelow, minPropPerGroup) clrs <- c("#E78AC3", "#8DA0CB") plotReceptLoss(exprMatrNml, exprMatrTum, rl, geneName="g7", clrs=clrs)
This function allows you to identify genes with loss of expression
receptLoss(exprMatrNml, exprMatrTum, nSdBelow, minPropPerGroup)
receptLoss(exprMatrNml, exprMatrTum, nSdBelow, minPropPerGroup)
exprMatrNml |
A matrix of expression values from normal tissue. Each row is a gene, and each column is a patient or sample. Genes should be in same order as exprMatrTum. |
exprMatrTum |
A matrix of expression values from tumor tissue. Each row is a gene, and each column is a patient or sample. Genes should be in same order as exprMatrNml. |
nSdBelow |
The number of SD below the mean of the adjacent normal tissue to set the boundary between tumor subgroups. |
minPropPerGroup |
A value between 0-1 that represents the minimum proportion of samples that should be present in each of the two subgroups (defined by the boundary set by nSdBelow) for a particular gene. |
a nx7 matrix, with n equaling the number of genes. The columns are as follows:
geneNm - the gene name
lowerBound - the lower bound, or the value 'nSdBelow' the mean of the normal tissue expression data.
propTumLessThBound - the proportion of tumor samples with expression levels less than 'lowerBound'
muAb - "mu above", the mean expression value of tumors greater than (ie above) the 'lowerBound'.
'muBl' - "mu below", the mean expression value of tumors less than (ie below) the 'lowerBound'.
'deltaMu' - the difference between 'muAb' and 'muBl'.
meetsMinPropPerGrp - a logical indicating whether the proportion of samples in each group is greater than that set by 'minPropPerGroup'.
exprMatrNml <- matrix(abs(rnorm(100, mean=2)), nrow=10) exprMatrTum <- matrix(abs(rnorm(100)), nrow=10) geneNames <- paste0(letters[seq_len(nrow(exprMatrNml))], seq_len(nrow(exprMatrNml))) rownames(exprMatrNml) <- rownames(exprMatrTum) <- geneNames nSdBelow <- 2 minPropPerGroup <- .2 rl <- receptLoss(exprMatrNml, exprMatrTum, nSdBelow, minPropPerGroup) head(rl)
exprMatrNml <- matrix(abs(rnorm(100, mean=2)), nrow=10) exprMatrTum <- matrix(abs(rnorm(100)), nrow=10) geneNames <- paste0(letters[seq_len(nrow(exprMatrNml))], seq_len(nrow(exprMatrNml))) rownames(exprMatrNml) <- rownames(exprMatrTum) <- geneNames nSdBelow <- 2 minPropPerGroup <- .2 rl <- receptLoss(exprMatrNml, exprMatrTum, nSdBelow, minPropPerGroup) head(rl)