Title: | Two-sample tests on a graph |
---|---|
Description: | DEGraph implements recent hypothesis testing methods which directly assess whether a particular gene network is differentially expressed between two conditions. This is to be contrasted with the more classical two-step approaches which first test individual genes, then test gene sets for enrichment in differentially expressed genes. These recent methods take into account the topology of the network to yield more powerful detection procedures. DEGraph provides methods to easily test all KEGG pathways for differential expression on any gene expression data set and tools to visualize the results. |
Authors: | Laurent Jacob, Pierre Neuvial and Sandrine Dudoit |
Maintainer: | Laurent Jacob <[email protected]> |
License: | GPL-3 |
Version: | 1.59.0 |
Built: | 2024-10-30 05:21:21 UTC |
Source: | https://github.com/bioc/DEGraph |
Performs the Adaptive Neyman test of Fan and Lin (1998).
AN.test(X1, X2, candK=1:ncol(X1), na.rm=FALSE)
AN.test(X1, X2, candK=1:ncol(X1), na.rm=FALSE)
X1 |
A n1 x p |
X2 |
A n2 x p |
candK |
A |
na.rm |
A |
A list
with class "htest" containing the following components:
A numeric
value, the test statistic.
A numeric
value, the corresponding p-value.
A numeric
value, the estimated true number of Fourier
components.
Laurent Jacob, Pierre Neuvial and Sandrine Dudoit
BS.test
()
graph.T2.test
()
hyper.test
()
library("KEGGgraph") ## library("NCIgraph") library("rrcov") data("Loi2008_DEGraphVignette") exprData <- exprLoi2008 classData <- classLoi2008 rn <- rownames(exprData) ## Retrieve expression levels data for genes from one KEGG pathway gr <- grListKEGG[[1]] gids <- translateKEGGID2GeneID(nodes(gr)) mm <- match(gids, rownames(exprData)) ## Keep genes from the graph that are present in the expression data set idxs <- which(!is.na(mm)) gr <- subGraph(nodes(gr)[idxs], gr) idxs <- which(is.na(mm)) if(length(idxs)) { print("Gene ID not found in expression data: ") str(gids[idxs]) } dat <- exprData[na.omit(mm), ] str(dat) X1 <- t(dat[, classData==0]) X2 <- t(dat[, classData==1]) ## DEGraph T2 test res <- testOneGraph(gr, exprData, classData, verbose=TRUE, prop=0.2) ## T2 test (Hotelling) rT2 <- T2.test(X1, X2) str(rT2) ## Adaptive Neyman test rAN <- AN.test(X1, X2, na.rm=TRUE) str(rAN) ## Adaptive Neyman test from Fan and Lin (1998) rAN <- AN.test(X1, X2, na.rm=TRUE) str(rAN) ## Test from Bai and Saranadasa (1996) rBS <- BS.test(X1, X2, na.rm=TRUE) str(rBS) ## Hypergeometric test pValues <- apply(exprData, 1, FUN=function(x) { tt <- t.test(x[classData==0], x[classData==1]) tt$p.value }) str(pValues) names(pValues) <- rownames(exprData) rHyper <- hyper.test(pValues, gids, thr=0.01) str(rHyper)
library("KEGGgraph") ## library("NCIgraph") library("rrcov") data("Loi2008_DEGraphVignette") exprData <- exprLoi2008 classData <- classLoi2008 rn <- rownames(exprData) ## Retrieve expression levels data for genes from one KEGG pathway gr <- grListKEGG[[1]] gids <- translateKEGGID2GeneID(nodes(gr)) mm <- match(gids, rownames(exprData)) ## Keep genes from the graph that are present in the expression data set idxs <- which(!is.na(mm)) gr <- subGraph(nodes(gr)[idxs], gr) idxs <- which(is.na(mm)) if(length(idxs)) { print("Gene ID not found in expression data: ") str(gids[idxs]) } dat <- exprData[na.omit(mm), ] str(dat) X1 <- t(dat[, classData==0]) X2 <- t(dat[, classData==1]) ## DEGraph T2 test res <- testOneGraph(gr, exprData, classData, verbose=TRUE, prop=0.2) ## T2 test (Hotelling) rT2 <- T2.test(X1, X2) str(rT2) ## Adaptive Neyman test rAN <- AN.test(X1, X2, na.rm=TRUE) str(rAN) ## Adaptive Neyman test from Fan and Lin (1998) rAN <- AN.test(X1, X2, na.rm=TRUE) str(rAN) ## Test from Bai and Saranadasa (1996) rBS <- BS.test(X1, X2, na.rm=TRUE) str(rBS) ## Hypergeometric test pValues <- apply(exprData, 1, FUN=function(x) { tt <- t.test(x[classData==0], x[classData==1]) tt$p.value }) str(pValues) names(pValues) <- rownames(exprData) rHyper <- hyper.test(pValues, gids, thr=0.01) str(rHyper)
This data set gives NCBI, Hugo and alternative gene symbols along with the cytoband and description for the 227 genes used in the DEGraph package vignette. This comes from the 15737 gene, 255 patient dataset of Loi et al. (2008) which was used to study resistance to tamoxifen treatment in hormone-dependent breast cancer.
annLoi2008
annLoi2008
A matrix of 227 lines and 5 columns.
Laurent Jacob, Pierre Neuvial and Sandrine Dudoit
Loi et al., Predicting prognosis using molecular profiling in estrogen receptor-positive breast cancer treated with tamoxifen. BMC Genomics, 9(1):239, 2008.
Loi et al., Predicting prognosis using molecular profiling in estrogen receptor-positive breast cancer treated with tamoxifen. BMC Genomics, 9(1):239, 2008.
data("Loi2008_DEGraphVignette") dim(annLoi2008) head(annLoi2008)
data("Loi2008_DEGraphVignette") dim(annLoi2008) head(annLoi2008)
Performs the test of Bai and Saranadasa (1996).
BS.test(X1, X2, na.rm=FALSE)
BS.test(X1, X2, na.rm=FALSE)
X1 |
A n1 x p |
X2 |
A n2 x p |
na.rm |
A |
A list
with class "htest" containing the following components:
Laurent Jacob, Pierre Neuvial and Sandrine Dudoit
AN.test
()
graph.T2.test
()
hyper.test
()
library("KEGGgraph") ## library("NCIgraph") library("rrcov") data("Loi2008_DEGraphVignette") exprData <- exprLoi2008 classData <- classLoi2008 rn <- rownames(exprData) ## Retrieve expression levels data for genes from one KEGG pathway gr <- grListKEGG[[1]] gids <- translateKEGGID2GeneID(nodes(gr)) mm <- match(gids, rownames(exprData)) ## Keep genes from the graph that are present in the expression data set idxs <- which(!is.na(mm)) gr <- subGraph(nodes(gr)[idxs], gr) idxs <- which(is.na(mm)) if(length(idxs)) { print("Gene ID not found in expression data: ") str(gids[idxs]) } dat <- exprData[na.omit(mm), ] str(dat) X1 <- t(dat[, classData==0]) X2 <- t(dat[, classData==1]) ## DEGraph T2 test res <- testOneGraph(gr, exprData, classData, verbose=TRUE, prop=0.2) ## T2 test (Hotelling) rT2 <- T2.test(X1, X2) str(rT2) ## Adaptive Neyman test rAN <- AN.test(X1, X2, na.rm=TRUE) str(rAN) ## Adaptive Neyman test from Fan and Lin (1998) rAN <- AN.test(X1, X2, na.rm=TRUE) str(rAN) ## Test from Bai and Saranadasa (1996) rBS <- BS.test(X1, X2, na.rm=TRUE) str(rBS) ## Hypergeometric test pValues <- apply(exprData, 1, FUN=function(x) { tt <- t.test(x[classData==0], x[classData==1]) tt$p.value }) str(pValues) names(pValues) <- rownames(exprData) rHyper <- hyper.test(pValues, gids, thr=0.01) str(rHyper)
library("KEGGgraph") ## library("NCIgraph") library("rrcov") data("Loi2008_DEGraphVignette") exprData <- exprLoi2008 classData <- classLoi2008 rn <- rownames(exprData) ## Retrieve expression levels data for genes from one KEGG pathway gr <- grListKEGG[[1]] gids <- translateKEGGID2GeneID(nodes(gr)) mm <- match(gids, rownames(exprData)) ## Keep genes from the graph that are present in the expression data set idxs <- which(!is.na(mm)) gr <- subGraph(nodes(gr)[idxs], gr) idxs <- which(is.na(mm)) if(length(idxs)) { print("Gene ID not found in expression data: ") str(gids[idxs]) } dat <- exprData[na.omit(mm), ] str(dat) X1 <- t(dat[, classData==0]) X2 <- t(dat[, classData==1]) ## DEGraph T2 test res <- testOneGraph(gr, exprData, classData, verbose=TRUE, prop=0.2) ## T2 test (Hotelling) rT2 <- T2.test(X1, X2) str(rT2) ## Adaptive Neyman test rAN <- AN.test(X1, X2, na.rm=TRUE) str(rAN) ## Adaptive Neyman test from Fan and Lin (1998) rAN <- AN.test(X1, X2, na.rm=TRUE) str(rAN) ## Test from Bai and Saranadasa (1996) rBS <- BS.test(X1, X2, na.rm=TRUE) str(rBS) ## Hypergeometric test pValues <- apply(exprData, 1, FUN=function(x) { tt <- t.test(x[classData==0], x[classData==1]) tt$p.value }) str(pValues) names(pValues) <- rownames(exprData) rHyper <- hyper.test(pValues, gids, thr=0.01) str(rHyper)
This data set gives resistance status data for the 255 patients used in the DEGraph package vignette. This comes from the 15737 gene, 255 patient dataset of Loi et al. (2008) which was used to study resistance to tamoxifen treatment in hormone-dependent breast cancer.
classLoi2008
classLoi2008
A vector of 255 elements which are either 0 (resistance to treatment) or 1 (sensitivity to treatment).
Laurent Jacob, Pierre Neuvial and Sandrine Dudoit
Loi et al., Predicting prognosis using molecular profiling in estrogen receptor-positive breast cancer treated with tamoxifen. BMC Genomics, 9(1):239, 2008.
Loi et al., Predicting prognosis using molecular profiling in estrogen receptor-positive breast cancer treated with tamoxifen. BMC Genomics, 9(1):239, 2008.
data("Loi2008_DEGraphVignette") dim(classLoi2008) head(classLoi2008)
data("Loi2008_DEGraphVignette") dim(classLoi2008) head(classLoi2008)
This data set gives gene expression data for a subset of 227 genes used in the DEGraph package vignette. This comes from the 15737 gene, 255 patient dataset of Loi et al. (2008) which was used to study resistance to tamoxifen treatment in hormone-dependent breast cancer.
exprLoi2008
exprLoi2008
A matrix of 227 lines and 255 columns.
The original data set corresponds to data processed by RMA and median-centered as available from the GSE6532 GEO archive: http://www.ncbi.nlm.nih.gov/geo/query/acc.cgi?acc=GSE6532.
These data were summarized from the probe set level to the gene level as follows. The expression level of a gene was defined as the expression level of the probe set with largest alignment score among all probe sets mapping to this gene according to the annotation in GSE6532. When the largest alignment score was achieved by several probe sets, the median expression level of those probe sets was taken.
Laurent Jacob, Pierre Neuvial and Sandrine Dudoit
Loi et al., Predicting prognosis using molecular profiling in estrogen receptor-positive breast cancer treated with tamoxifen. BMC Genomics, 9(1):239, 2008.
Loi et al., Predicting prognosis using molecular profiling in estrogen receptor-positive breast cancer treated with tamoxifen. BMC Genomics, 9(1):239, 2008.
data("Loi2008_DEGraphVignette") dim(exprLoi2008) head(exprLoi2008)
data("Loi2008_DEGraphVignette") dim(exprLoi2008) head(exprLoi2008)
Given a graph, returns a list of its connected components (which are also graph objects), ordered by decreasing number of nodes.
getConnectedComponentList(graph, verbose=FALSE)
getConnectedComponentList(graph, verbose=FALSE)
graph |
A |
verbose |
If |
A list
containing a graph
object for each connected
component of the input graph, ordered by decreasing number of nodes
Laurent Jacob, Pierre Neuvial and Sandrine Dudoit
data("Loi2008_DEGraphVignette") exprData <- exprLoi2008 rn <- rownames(exprData) ## Retrieve expression levels data for genes from one KEGG pathway graph <- grListKEGG[[1]] pname <- attr(graph, "label") cat(verbose, "Pathway name: ", pname) sgraph <- getSignedGraph(graph, verbose=TRUE) print(sgraph) graphList <- getConnectedComponentList(graph, verbose=TRUE) print(graphList)
data("Loi2008_DEGraphVignette") exprData <- exprLoi2008 rn <- rownames(exprData) ## Retrieve expression levels data for genes from one KEGG pathway graph <- grListKEGG[[1]] pname <- attr(graph, "label") cat(verbose, "Pathway name: ", pname) sgraph <- getSignedGraph(graph, verbose=TRUE) print(sgraph) graphList <- getConnectedComponentList(graph, verbose=TRUE) print(graphList)
Builds a graph for each of the KEGG pathways.
getKEGGPathways(path=NULL, rootPath="networkData/ftp.genome.jp/pub/kegg/xml/kgml", organism="hsa", metaTag=c("non-metabolic", "metabolic"), pattern=NULL, verbose=FALSE)
getKEGGPathways(path=NULL, rootPath="networkData/ftp.genome.jp/pub/kegg/xml/kgml", organism="hsa", metaTag=c("non-metabolic", "metabolic"), pattern=NULL, verbose=FALSE)
path |
A |
rootPath |
A |
organism |
A |
metaTag |
A |
pattern |
An optional |
verbose |
If |
If 'path' is supplied, KGML files in this directory are loaded. Otherwise, KGML files are assumed to be in <rootPath>/<metaTag>/"organisms"/<organism>, which mirrors the structure of the KEGG KGML file repository.
A list
containing a graph
object for each KEGG pathway with at least one edge.
Laurent Jacob, Pierre Neuvial and Sandrine Dudoit
library("Rgraphviz") library("KEGGgraph") ## example of KGML files path <- system.file("extdata", package="KEGGgraph") grList <- getKEGGPathways(path=path, verbose=TRUE) print(grList) graph <- grList[[1]] plotKEGGgraph(graph) ## Not run: ## Download all human KGML pathways locally pathname <- system.file("downloadScripts", "downloadKeggXmlFiles.R", package="DEGraph") source(pathname) ## Load some of them grList <- getKEGGPathways(pattern="040", verbose=TRUE) print(grList) graph <- grList[[1]] plotKEGGgraph(graph) ## End(Not run)
library("Rgraphviz") library("KEGGgraph") ## example of KGML files path <- system.file("extdata", package="KEGGgraph") grList <- getKEGGPathways(path=path, verbose=TRUE) print(grList) graph <- grList[[1]] plotKEGGgraph(graph) ## Not run: ## Download all human KGML pathways locally pathname <- system.file("downloadScripts", "downloadKeggXmlFiles.R", package="DEGraph") source(pathname) ## Load some of them grList <- getKEGGPathways(pattern="040", verbose=TRUE) print(grList) graph <- grList[[1]] plotKEGGgraph(graph) ## End(Not run)
Given a graph, builds a signed version of the adjacency matrix taking into account the type of interaction (e.g., activation or inhibition).
getSignedGraph(graph, positiveInteractionLabels=c("activation", "expression"), negativeInteractionLabels=c("inhibition", "repression"), verbose=FALSE)
getSignedGraph(graph, positiveInteractionLabels=c("activation", "expression"), negativeInteractionLabels=c("inhibition", "repression"), verbose=FALSE)
graph |
A |
positiveInteractionLabels |
A |
negativeInteractionLabels |
A |
verbose |
If |
This function returns a squared matrix whose (i,j) entry is:
if edges i and j are not connected
if edges i and j are connected by a positive interaction
if edges i and j are connected by a negative interaction.
By construction, the absolute value of this matrix is the adjacency matrix of the graph. Edges which cannot interpreted as corresponding to a positive or a negative interaction are marked as not connected.
Laurent Jacob, Pierre Neuvial and Sandrine Dudoit
data("Loi2008_DEGraphVignette") exprData <- exprLoi2008 rn <- rownames(exprData) ## Retrieve expression levels data for genes from one KEGG pathway graph <- grListKEGG[[1]] pname <- attr(graph, "label") cat(verbose, "Pathway name: ", pname) sgraph <- getSignedGraph(graph, verbose=TRUE) print(sgraph) graphList <- getConnectedComponentList(graph, verbose=TRUE) print(graphList)
data("Loi2008_DEGraphVignette") exprData <- exprLoi2008 rn <- rownames(exprData) ## Retrieve expression levels data for genes from one KEGG pathway graph <- grListKEGG[[1]] pname <- attr(graph, "label") cat(verbose, "Pathway name: ", pname) sgraph <- getSignedGraph(graph, verbose=TRUE) print(sgraph) graphList <- getConnectedComponentList(graph, verbose=TRUE) print(graphList)
Performs the Hotelling T2 test in Fourier space.
graph.T2.test(X1, X2, G=NULL, lfA=NULL, ..., k=ncol(X1))
graph.T2.test(X1, X2, G=NULL, lfA=NULL, ..., k=ncol(X1))
X1 |
A n1 x p |
X2 |
A n2 x p |
G |
An object of class |
lfA |
A list returned by |
... |
Further arguments to be passed to |
k |
A |
A list
with class "htest", as returned by T2.test
.
Laurent Jacob, Pierre Neuvial and Sandrine Dudoit
library("rrcov") ## Some parameters n1 <- n2 <- 20 nnodes <- nedges <- 20 k <- 3 ncp <- 0.5 sigma <- diag(nnodes)/sqrt(nnodes) ## Build graph, decompose laplacian G <- randomWAMGraph(nnodes=nnodes,nedges=nedges) A <- G@adjMat lfA <- laplacianFromA(A,ltype="unnormalized") U <- lfA$U l <- lfA$l ## Build two samples with smooth mean shift X <- twoSampleFromGraph(n1,n2,shiftM2=ncp,sigma,U=U,k=k) ## Do hypothesis testing t <- T2.test(X$X1,X$X2) # Raw T-square print(t$p.value) tu <- graph.T2.test(X$X1,X$X2,lfA=lfA,k=k) # Filtered T-squares print(tu$p.value)
library("rrcov") ## Some parameters n1 <- n2 <- 20 nnodes <- nedges <- 20 k <- 3 ncp <- 0.5 sigma <- diag(nnodes)/sqrt(nnodes) ## Build graph, decompose laplacian G <- randomWAMGraph(nnodes=nnodes,nedges=nedges) A <- G@adjMat lfA <- laplacianFromA(A,ltype="unnormalized") U <- lfA$U l <- lfA$l ## Build two samples with smooth mean shift X <- twoSampleFromGraph(n1,n2,shiftM2=ncp,sigma,U=U,k=k) ## Do hypothesis testing t <- T2.test(X$X1,X$X2) # Raw T-square print(t$p.value) tu <- graph.T2.test(X$X1,X$X2,lfA=lfA,k=k) # Filtered T-squares print(tu$p.value)
This data set gives KEGGgraph objects for two KEGG non-metabolic pathways ("Natural killer cell mediated cytotoxicity" and "Insulin signaling pathway").
grListKEGG
grListKEGG
A list of two elements.
Laurent Jacob, Pierre Neuvial and Sandrine Dudoit
library("Rgraphviz") data("Loi2008_DEGraphVignette") grListKEGG plot(grListKEGG[[1]])
library("Rgraphviz") data("Loi2008_DEGraphVignette") grListKEGG plot(grListKEGG[[1]])
Performs an hypergeometric test of enrichment of a set of hypotheses in significant elements.
hyper.test(p.values, testSet, thr=0.001, universe=length(p.values), verbose=FALSE)
hyper.test(p.values, testSet, thr=0.001, universe=length(p.values), verbose=FALSE)
p.values |
A named |
testSet |
A |
thr |
A |
universe |
An |
verbose |
If |
A list
with class "htest" containing the following components:
Laurent Jacob, Pierre Neuvial and Sandrine Dudoit
AN.test
()
BS.test
()
graph.T2.test
()
library("KEGGgraph") ## library("NCIgraph") library("rrcov") data("Loi2008_DEGraphVignette") exprData <- exprLoi2008 classData <- classLoi2008 rn <- rownames(exprData) ## Retrieve expression levels data for genes from one KEGG pathway gr <- grListKEGG[[1]] gids <- translateKEGGID2GeneID(nodes(gr)) mm <- match(gids, rownames(exprData)) ## Keep genes from the graph that are present in the expression data set idxs <- which(!is.na(mm)) gr <- subGraph(nodes(gr)[idxs], gr) idxs <- which(is.na(mm)) if(length(idxs)) { print("Gene ID not found in expression data: ") str(gids[idxs]) } dat <- exprData[na.omit(mm), ] str(dat) X1 <- t(dat[, classData==0]) X2 <- t(dat[, classData==1]) ## DEGraph T2 test res <- testOneGraph(gr, exprData, classData, verbose=TRUE, prop=0.2) ## T2 test (Hotelling) rT2 <- T2.test(X1, X2) str(rT2) ## Adaptive Neyman test rAN <- AN.test(X1, X2, na.rm=TRUE) str(rAN) ## Adaptive Neyman test from Fan and Lin (1998) rAN <- AN.test(X1, X2, na.rm=TRUE) str(rAN) ## Test from Bai and Saranadasa (1996) rBS <- BS.test(X1, X2, na.rm=TRUE) str(rBS) ## Hypergeometric test pValues <- apply(exprData, 1, FUN=function(x) { tt <- t.test(x[classData==0], x[classData==1]) tt$p.value }) str(pValues) names(pValues) <- rownames(exprData) rHyper <- hyper.test(pValues, gids, thr=0.01) str(rHyper)
library("KEGGgraph") ## library("NCIgraph") library("rrcov") data("Loi2008_DEGraphVignette") exprData <- exprLoi2008 classData <- classLoi2008 rn <- rownames(exprData) ## Retrieve expression levels data for genes from one KEGG pathway gr <- grListKEGG[[1]] gids <- translateKEGGID2GeneID(nodes(gr)) mm <- match(gids, rownames(exprData)) ## Keep genes from the graph that are present in the expression data set idxs <- which(!is.na(mm)) gr <- subGraph(nodes(gr)[idxs], gr) idxs <- which(is.na(mm)) if(length(idxs)) { print("Gene ID not found in expression data: ") str(gids[idxs]) } dat <- exprData[na.omit(mm), ] str(dat) X1 <- t(dat[, classData==0]) X2 <- t(dat[, classData==1]) ## DEGraph T2 test res <- testOneGraph(gr, exprData, classData, verbose=TRUE, prop=0.2) ## T2 test (Hotelling) rT2 <- T2.test(X1, X2) str(rT2) ## Adaptive Neyman test rAN <- AN.test(X1, X2, na.rm=TRUE) str(rAN) ## Adaptive Neyman test from Fan and Lin (1998) rAN <- AN.test(X1, X2, na.rm=TRUE) str(rAN) ## Test from Bai and Saranadasa (1996) rBS <- BS.test(X1, X2, na.rm=TRUE) str(rBS) ## Hypergeometric test pValues <- apply(exprData, 1, FUN=function(x) { tt <- t.test(x[classData==0], x[classData==1]) tt$p.value }) str(pValues) names(pValues) <- rownames(exprData) rHyper <- hyper.test(pValues, gids, thr=0.01) str(rHyper)
Calculates the Laplacian associated to an adjacency matrix.
laplacianFromA(A, k=1, ltype=c("meanInfluence", "normalized", "unnormalized", "totalInfluence"))
laplacianFromA(A, k=1, ltype=c("meanInfluence", "normalized", "unnormalized", "totalInfluence"))
A |
The adjacency matrix of the graph. |
k |
... |
ltype |
A |
A list
containing the following components:
Eigenvectors of the graph Laplacian.
Eigenvalues of the graph Laplacian
Multiplicity of '0' as eigenvalue.
Laurent Jacob, Pierre Neuvial and Sandrine Dudoit
library("KEGGgraph") library("rrcov") ## Create a random graph graph <- randomWAMGraph(nnodes=5, nedges=7, verbose=TRUE) plot(graph) ## Retrieve its adjacency matrix A <- graph@adjMat ## write it to KGML file grPathname <- "randomWAMGraph.xml" writeAdjacencyMatrix2KGML(A, pathname=grPathname, verbose=TRUE, overwrite=TRUE) ## read it from file gr <- parseKGML2Graph(grPathname) ## Two examples of Laplacians from the same graph lapMI <- laplacianFromA(A, ltype="meanInfluence") print(lapMI) lapN <- laplacianFromA(A, ltype="normalized") print(lapN) U <- lapN$U p <- nrow(A) sigma <- diag(p)/sqrt(p) X <- twoSampleFromGraph(100, 120, shiftM2=1, sigma, U=U, k=3) ## T2 t <- T2.test(X$X1,X$X2) str(t) tu <- graph.T2.test(X$X1, X$X2, lfA=lapMI, k=3) str(tu)
library("KEGGgraph") library("rrcov") ## Create a random graph graph <- randomWAMGraph(nnodes=5, nedges=7, verbose=TRUE) plot(graph) ## Retrieve its adjacency matrix A <- graph@adjMat ## write it to KGML file grPathname <- "randomWAMGraph.xml" writeAdjacencyMatrix2KGML(A, pathname=grPathname, verbose=TRUE, overwrite=TRUE) ## read it from file gr <- parseKGML2Graph(grPathname) ## Two examples of Laplacians from the same graph lapMI <- laplacianFromA(A, ltype="meanInfluence") print(lapMI) lapN <- laplacianFromA(A, ltype="normalized") print(lapN) U <- lapN$U p <- nrow(A) sigma <- diag(p)/sqrt(p) X <- twoSampleFromGraph(100, 120, shiftM2=1, sigma, U=U, k=3) ## T2 t <- T2.test(X$X1,X$X2) str(t) tu <- graph.T2.test(X$X1, X$X2, lfA=lapMI, k=3) str(tu)
Plots a graph with nodes colored according to a quantitative variable.
plotValuedGraph(graph, values=NULL, nodeLabels=nodes(graph), qMax=0.95, colorPalette=heat.colors(10), adjustColorRange=FALSE, symmetrizeArrows=FALSE, height=1, lwd=1, cex=1, ..., verbose=FALSE)
plotValuedGraph(graph, values=NULL, nodeLabels=nodes(graph), qMax=0.95, colorPalette=heat.colors(10), adjustColorRange=FALSE, symmetrizeArrows=FALSE, height=1, lwd=1, cex=1, ..., verbose=FALSE)
graph |
A |
values |
A named |
nodeLabels |
A |
qMax |
A |
colorPalette |
A |
adjustColorRange |
A |
symmetrizeArrows |
A |
height |
A |
lwd |
A |
cex |
A |
... |
Further arguments to be passed to 'edgeRenderInfo' and 'nodeRenderInfo'. |
verbose |
If |
A list
containing the following components:
The 'graph' object as plotted.
The break points in the supplied values (can be used for plotting a legend).
Laurent Jacob, Pierre Neuvial and Sandrine Dudoit
library("Rgraphviz") library("KEGGgraph") ## library("NCIgraph") data("Loi2008_DEGraphVignette") exprData <- exprLoi2008 classData <- classLoi2008 annData <- annLoi2008 rn <- rownames(exprData) ## Retrieve expression levels data for genes from one KEGG pathway graph <- grListKEGG[[1]] pname <- attr(graph, "label") print(pname) ## DEGraph T2 test resList <- testOneGraph(graph, exprData, classData, verbose=TRUE, prop=0.2) ## Largest connected component res <- resList[[1]] gr <- res$graph ## individual t statistics shift <- apply(exprData, 1, FUN=function(x) { tt <- t.test(x[classData==0], x[classData==1]) tt$statistic }) names(shift) <- translateGeneID2KEGGID(names(shift)) ## color palette if (require(marray)) { pal <- maPalette(low="red", high="green", mid="black", k=100) } else { pal <- heat.colors(100) } ## plot results dn <- getDisplayName(gr, shortLabel=TRUE) mm <- match(translateKEGGID2GeneID(nodes(gr)), rownames(annData)) dn <- annData[mm, "NCBI.gene.symbol"] pvg <- plotValuedGraph(gr, values=shift, nodeLabels=dn, qMax=0.95, colorPalette=pal, height=40, lwd=1, verbose=TRUE, cex=0.5) title(pname) txt1 <- sprintf("p(T2)=%s", signif(res$p.value[1], 2)) txt2 <- sprintf("p(T2F[%s])=%s", res$k, signif(res$p.value[2])) txt <- paste(txt1, txt2, sep="\n") stext(side=3, pos=1, txt) if (require(fields)) { image.plot(legend.only=TRUE, zlim=range(pvg$breaks), col=pal, legend.shrink=0.3, legend.width=0.8, legend.lab="t-scores", legend.mar=3.3) }
library("Rgraphviz") library("KEGGgraph") ## library("NCIgraph") data("Loi2008_DEGraphVignette") exprData <- exprLoi2008 classData <- classLoi2008 annData <- annLoi2008 rn <- rownames(exprData) ## Retrieve expression levels data for genes from one KEGG pathway graph <- grListKEGG[[1]] pname <- attr(graph, "label") print(pname) ## DEGraph T2 test resList <- testOneGraph(graph, exprData, classData, verbose=TRUE, prop=0.2) ## Largest connected component res <- resList[[1]] gr <- res$graph ## individual t statistics shift <- apply(exprData, 1, FUN=function(x) { tt <- t.test(x[classData==0], x[classData==1]) tt$statistic }) names(shift) <- translateGeneID2KEGGID(names(shift)) ## color palette if (require(marray)) { pal <- maPalette(low="red", high="green", mid="black", k=100) } else { pal <- heat.colors(100) } ## plot results dn <- getDisplayName(gr, shortLabel=TRUE) mm <- match(translateKEGGID2GeneID(nodes(gr)), rownames(annData)) dn <- annData[mm, "NCBI.gene.symbol"] pvg <- plotValuedGraph(gr, values=shift, nodeLabels=dn, qMax=0.95, colorPalette=pal, height=40, lwd=1, verbose=TRUE, cex=0.5) title(pname) txt1 <- sprintf("p(T2)=%s", signif(res$p.value[1], 2)) txt2 <- sprintf("p(T2F[%s])=%s", res$k, signif(res$p.value[2])) txt <- paste(txt1, txt2, sep="\n") stext(side=3, pos=1, txt) if (require(fields)) { image.plot(legend.only=TRUE, zlim=range(pvg$breaks), col=pal, legend.shrink=0.3, legend.width=0.8, legend.lab="t-scores", legend.mar=3.3) }
Generates a random graph.
randomWAMGraph(nnodes=5, nedges=nnodes, verbose=FALSE)
randomWAMGraph(nnodes=5, nedges=nnodes, verbose=FALSE)
nnodes |
A |
nedges |
A |
verbose |
If |
An object of class graphAM
.
Laurent Jacob, Pierre Neuvial and Sandrine Dudoit
library("KEGGgraph") library("rrcov") ## Create a random graph graph <- randomWAMGraph(nnodes=5, nedges=7, verbose=TRUE) plot(graph) ## Retrieve its adjacency matrix A <- graph@adjMat ## write it to KGML file grPathname <- "randomWAMGraph.xml" writeAdjacencyMatrix2KGML(A, pathname=grPathname, verbose=TRUE, overwrite=TRUE) ## read it from file gr <- parseKGML2Graph(grPathname) ## Two examples of Laplacians from the same graph lapMI <- laplacianFromA(A, ltype="meanInfluence") print(lapMI) lapN <- laplacianFromA(A, ltype="normalized") print(lapN) U <- lapN$U p <- nrow(A) sigma <- diag(p)/sqrt(p) X <- twoSampleFromGraph(100, 120, shiftM2=1, sigma, U=U, k=3) ## T2 t <- T2.test(X$X1,X$X2) str(t) tu <- graph.T2.test(X$X1, X$X2, lfA=lapMI, k=3) str(tu)
library("KEGGgraph") library("rrcov") ## Create a random graph graph <- randomWAMGraph(nnodes=5, nedges=7, verbose=TRUE) plot(graph) ## Retrieve its adjacency matrix A <- graph@adjMat ## write it to KGML file grPathname <- "randomWAMGraph.xml" writeAdjacencyMatrix2KGML(A, pathname=grPathname, verbose=TRUE, overwrite=TRUE) ## read it from file gr <- parseKGML2Graph(grPathname) ## Two examples of Laplacians from the same graph lapMI <- laplacianFromA(A, ltype="meanInfluence") print(lapMI) lapN <- laplacianFromA(A, ltype="normalized") print(lapN) U <- lapN$U p <- nrow(A) sigma <- diag(p)/sqrt(p) X <- twoSampleFromGraph(100, 120, shiftM2=1, sigma, U=U, k=3) ## T2 t <- T2.test(X$X1,X$X2) str(t) tu <- graph.T2.test(X$X1, X$X2, lfA=lapMI, k=3) str(tu)
Applies a series of two-sample tests to a connected graph using various statistics.
testOneConnectedComponent(graph, data, classes, ..., prop=0.2, verbose=FALSE)
testOneConnectedComponent(graph, data, classes, ..., prop=0.2, verbose=FALSE)
graph |
A |
data |
A ' |
classes |
|
... |
Further arguments to be passed to |
prop |
A |
verbose |
If |
This function performs the test, assuming that all genes in the graph are represented in the expression data set, in order not to have to modify the graph topology.
Interaction signs are used if available in the graph ('getSignedGraph' is not called here, in order not to have to modify the graph topology.).
The graph given as input has to have only one
connex component. It can be retrieved from the output of
getConnectedComponentList
().
A structured list
containing the p-values of the tests, the
graph
object of the connected component and the number
of retained Fourier dimensions.
Laurent Jacob, Pierre Neuvial and Sandrine Dudoit
testOneGraph
()
getConnectedComponentList
()
library("rrcov") ## Some parameters n1 <- n2 <- 20 nnodes <- nedges <- 20 k <- 3 ncp <- 0.5 sigma <- diag(nnodes)/sqrt(nnodes) ## Build graph, decompose laplacian G <- randomWAMGraph(nnodes=nnodes,nedges=nedges) A <- G@adjMat lfA <- laplacianFromA(A,ltype="unnormalized") U <- lfA$U l <- lfA$l ## Build two samples with smooth mean shift X <- twoSampleFromGraph(n1,n2,shiftM2=ncp,sigma,U=U,k=k) ## Do hypothesis testing t <- T2.test(X$X1,X$X2) # Raw T-square print(t$p.value) tu <- graph.T2.test(X$X1,X$X2,lfA=lfA,k=k) # Filtered T-squares print(tu$p.value)
library("rrcov") ## Some parameters n1 <- n2 <- 20 nnodes <- nedges <- 20 k <- 3 ncp <- 0.5 sigma <- diag(nnodes)/sqrt(nnodes) ## Build graph, decompose laplacian G <- randomWAMGraph(nnodes=nnodes,nedges=nedges) A <- G@adjMat lfA <- laplacianFromA(A,ltype="unnormalized") U <- lfA$U l <- lfA$l ## Build two samples with smooth mean shift X <- twoSampleFromGraph(n1,n2,shiftM2=ncp,sigma,U=U,k=k) ## Do hypothesis testing t <- T2.test(X$X1,X$X2) # Raw T-square print(t$p.value) tu <- graph.T2.test(X$X1,X$X2,lfA=lfA,k=k) # Filtered T-squares print(tu$p.value)
Applies a serie of two-sample tests to each connected component of a graph using various statistics.
testOneGraph(graph, data, classes, useInteractionSigns=TRUE, ..., verbose=FALSE)
testOneGraph(graph, data, classes, useInteractionSigns=TRUE, ..., verbose=FALSE)
graph |
A |
data |
A 'matrix' (size: number 'p' of genes x number 'n' of samples) of gene expression. |
classes |
A 'vector' (length: 'n') of class assignments. |
useInteractionSigns |
A |
... |
Further arguments to be passed to testOneConnectedComponent. |
verbose |
If |
A structured list
containing the p-values of the tests, the
graph
object of the connected component and the number of
retained Fourier dimensions.
Laurent Jacob, Pierre Neuvial and Sandrine Dudoit
library("Rgraphviz") library("KEGGgraph") ## library("NCIgraph") data("Loi2008_DEGraphVignette") exprData <- exprLoi2008 classData <- classLoi2008 annData <- annLoi2008 rn <- rownames(exprData) ## Retrieve expression levels data for genes from one KEGG pathway graph <- grListKEGG[[1]] pname <- attr(graph, "label") print(pname) ## DEGraph T2 test resList <- testOneGraph(graph, exprData, classData, verbose=TRUE, prop=0.2) ## Largest connected component res <- resList[[1]] gr <- res$graph ## individual t statistics shift <- apply(exprData, 1, FUN=function(x) { tt <- t.test(x[classData==0], x[classData==1]) tt$statistic }) names(shift) <- translateGeneID2KEGGID(names(shift)) ## color palette if (require(marray)) { pal <- maPalette(low="red", high="green", mid="black", k=100) } else { pal <- heat.colors(100) } ## plot results dn <- getDisplayName(gr, shortLabel=TRUE) mm <- match(translateKEGGID2GeneID(nodes(gr)), rownames(annData)) dn <- annData[mm, "NCBI.gene.symbol"] pvg <- plotValuedGraph(gr, values=shift, nodeLabels=dn, qMax=0.95, colorPalette=pal, height=40, lwd=1, verbose=TRUE, cex=0.5) title(pname) txt1 <- sprintf("p(T2)=%s", signif(res$p.value[1], 2)) txt2 <- sprintf("p(T2F[%s])=%s", res$k, signif(res$p.value[2])) txt <- paste(txt1, txt2, sep="\n") stext(side=3, pos=1, txt) if (require(fields)) { image.plot(legend.only=TRUE, zlim=range(pvg$breaks), col=pal, legend.shrink=0.3, legend.width=0.8, legend.lab="t-scores", legend.mar=3.3) }
library("Rgraphviz") library("KEGGgraph") ## library("NCIgraph") data("Loi2008_DEGraphVignette") exprData <- exprLoi2008 classData <- classLoi2008 annData <- annLoi2008 rn <- rownames(exprData) ## Retrieve expression levels data for genes from one KEGG pathway graph <- grListKEGG[[1]] pname <- attr(graph, "label") print(pname) ## DEGraph T2 test resList <- testOneGraph(graph, exprData, classData, verbose=TRUE, prop=0.2) ## Largest connected component res <- resList[[1]] gr <- res$graph ## individual t statistics shift <- apply(exprData, 1, FUN=function(x) { tt <- t.test(x[classData==0], x[classData==1]) tt$statistic }) names(shift) <- translateGeneID2KEGGID(names(shift)) ## color palette if (require(marray)) { pal <- maPalette(low="red", high="green", mid="black", k=100) } else { pal <- heat.colors(100) } ## plot results dn <- getDisplayName(gr, shortLabel=TRUE) mm <- match(translateKEGGID2GeneID(nodes(gr)), rownames(annData)) dn <- annData[mm, "NCBI.gene.symbol"] pvg <- plotValuedGraph(gr, values=shift, nodeLabels=dn, qMax=0.95, colorPalette=pal, height=40, lwd=1, verbose=TRUE, cex=0.5) title(pname) txt1 <- sprintf("p(T2)=%s", signif(res$p.value[1], 2)) txt2 <- sprintf("p(T2F[%s])=%s", res$k, signif(res$p.value[2])) txt <- paste(txt1, txt2, sep="\n") stext(side=3, pos=1, txt) if (require(fields)) { image.plot(legend.only=TRUE, zlim=range(pvg$breaks), col=pal, legend.shrink=0.3, legend.width=0.8, legend.lab="t-scores", legend.mar=3.3) }
Given a basis (typically the eigenvectors of a graph Laplacian), builds two multivariate normal samples with mean shift located in the first elements of the basis.
twoSampleFromGraph(n1=20, n2=n1, shiftM2=0, sigma, U, k=ceiling(ncol(U)/3))
twoSampleFromGraph(n1=20, n2=n1, shiftM2=0, sigma, U, k=ceiling(ncol(U)/3))
n1 |
An |
n2 |
An |
shiftM2 |
A |
sigma |
A matrix giving the covariance structure of each sample. |
U |
A matrix giving the desired basis. |
k |
An |
A list
with named elements:
The first sample in the original basis (before transformation by U).
The second sample in the original basis (before transformation by U).
The first sample in the specified basis (after transformation by U).
The second sample in the specified basis (after transformation by U).
The population mean of F1
The population mean of F2
mu1 - mu2
Laurent Jacob, Pierre Neuvial and Sandrine Dudoit
library("KEGGgraph") library("rrcov") ## Create a random graph graph <- randomWAMGraph(nnodes=5, nedges=7, verbose=TRUE) plot(graph) ## Retrieve its adjacency matrix A <- graph@adjMat ## write it to KGML file grPathname <- "randomWAMGraph.xml" writeAdjacencyMatrix2KGML(A, pathname=grPathname, verbose=TRUE, overwrite=TRUE) ## read it from file gr <- parseKGML2Graph(grPathname) ## Two examples of Laplacians from the same graph lapMI <- laplacianFromA(A, ltype="meanInfluence") print(lapMI) lapN <- laplacianFromA(A, ltype="normalized") print(lapN) U <- lapN$U p <- nrow(A) sigma <- diag(p)/sqrt(p) X <- twoSampleFromGraph(100, 120, shiftM2=1, sigma, U=U, k=3) ## T2 t <- T2.test(X$X1,X$X2) str(t) tu <- graph.T2.test(X$X1, X$X2, lfA=lapMI, k=3) str(tu)
library("KEGGgraph") library("rrcov") ## Create a random graph graph <- randomWAMGraph(nnodes=5, nedges=7, verbose=TRUE) plot(graph) ## Retrieve its adjacency matrix A <- graph@adjMat ## write it to KGML file grPathname <- "randomWAMGraph.xml" writeAdjacencyMatrix2KGML(A, pathname=grPathname, verbose=TRUE, overwrite=TRUE) ## read it from file gr <- parseKGML2Graph(grPathname) ## Two examples of Laplacians from the same graph lapMI <- laplacianFromA(A, ltype="meanInfluence") print(lapMI) lapN <- laplacianFromA(A, ltype="normalized") print(lapN) U <- lapN$U p <- nrow(A) sigma <- diag(p)/sqrt(p) X <- twoSampleFromGraph(100, 120, shiftM2=1, sigma, U=U, k=3) ## T2 t <- T2.test(X$X1,X$X2) str(t) tu <- graph.T2.test(X$X1, X$X2, lfA=lapMI, k=3) str(tu)
Writes an adjacency matrix into an XML file.
writeAdjacencyMatrix2KGML(mat, pathname, nodePrefix="n", overwrite=FALSE, ..., verbose=FALSE)
writeAdjacencyMatrix2KGML(mat, pathname, nodePrefix="n", overwrite=FALSE, ..., verbose=FALSE)
mat |
A |
pathname |
The full path name of the XML file to be written. |
nodePrefix |
A |
overwrite |
If |
... |
Further arguments to be passed to plotKEGGgraph. |
verbose |
If |
None.
Laurent Jacob, Pierre Neuvial and Sandrine Dudoit
library("KEGGgraph") library("rrcov") ## Create a random graph graph <- randomWAMGraph(nnodes=5, nedges=7, verbose=TRUE) plot(graph) ## Retrieve its adjacency matrix A <- graph@adjMat ## write it to KGML file grPathname <- "randomWAMGraph.xml" writeAdjacencyMatrix2KGML(A, pathname=grPathname, verbose=TRUE, overwrite=TRUE) ## read it from file gr <- parseKGML2Graph(grPathname) ## Two examples of Laplacians from the same graph lapMI <- laplacianFromA(A, ltype="meanInfluence") print(lapMI) lapN <- laplacianFromA(A, ltype="normalized") print(lapN) U <- lapN$U p <- nrow(A) sigma <- diag(p)/sqrt(p) X <- twoSampleFromGraph(100, 120, shiftM2=1, sigma, U=U, k=3) ## T2 t <- T2.test(X$X1,X$X2) str(t) tu <- graph.T2.test(X$X1, X$X2, lfA=lapMI, k=3) str(tu)
library("KEGGgraph") library("rrcov") ## Create a random graph graph <- randomWAMGraph(nnodes=5, nedges=7, verbose=TRUE) plot(graph) ## Retrieve its adjacency matrix A <- graph@adjMat ## write it to KGML file grPathname <- "randomWAMGraph.xml" writeAdjacencyMatrix2KGML(A, pathname=grPathname, verbose=TRUE, overwrite=TRUE) ## read it from file gr <- parseKGML2Graph(grPathname) ## Two examples of Laplacians from the same graph lapMI <- laplacianFromA(A, ltype="meanInfluence") print(lapMI) lapN <- laplacianFromA(A, ltype="normalized") print(lapN) U <- lapN$U p <- nrow(A) sigma <- diag(p)/sqrt(p) X <- twoSampleFromGraph(100, 120, shiftM2=1, sigma, U=U, k=3) ## T2 t <- T2.test(X$X1,X$X2) str(t) tu <- graph.T2.test(X$X1, X$X2, lfA=lapMI, k=3) str(tu)