Title: | Routines for the functional analysis of biological networks |
---|---|
Description: | This package provides functions for the integrated analysis of protein-protein interaction networks and the detection of functional modules. Different datasets can be integrated into the network by assigning p-values of statistical tests to the nodes of the network. E.g. p-values obtained from the differential expression of the genes from an Affymetrix array are assigned to the nodes of the network. By fitting a beta-uniform mixture model and calculating scores from the p-values, overall scores of network regions can be calculated and an integer linear programming algorithm identifies the maximum scoring subnetwork. |
Authors: | Marcus Dittrich and Daniela Beisser |
Maintainer: | Marcus Dittrich <[email protected]> |
License: | GPL (>= 2) |
Version: | 1.67.0 |
Built: | 2024-11-19 03:26:53 UTC |
Source: | https://github.com/bioc/BioNet |
This package provides functions for the integrated analysis of biological networks and the detection of functional modules. Different datasets can be integrated into the network by assigning p-values derived from statistical tests to the nodes of the network. E.g. p-values obtained from the differential expression of genes from an Affymetrix array are assigned to the nodes of an protein-protein interaction network. By fitting a beta-uniform mixture model and calculating scores from the p-values, overall scores of network regions can be calculated and an integer linear programming algorithm identifies the maximum scoring subnetwork.
Package: | BioNet |
Type: | Package |
Version: | 1.29.1 |
Date: | 2015-09-11 |
License: | GPL (>=2) |
LazyLoad: | yes |
Marcus Dittrich, Daniela Beisser
Maintainer: Marcus Dittrich <[email protected]>
M. T. Dittrich, G. W. Klau, A. Rosenwald, T. Dandekar and T. Mueller (2008) Identifying functional modules in protein-protein interaction networks: an integrated exact approach. (ISMB2008) Bioinformatics 24: 13. i223-i231 Jul.
D. Beisser, G. W. Klau, T. Dandekar, T. Mueller and M. Dittrich (2010) BioNet: an R-package for the Functional Analysis of Biological Networks. Bioinformatics 26:08. 1129-1130 Apr.
The function aggregates several p-values into one p-value of p-values based on the order statistics of p-values. An overall p-value is given by the ith order statistic.
aggrPvals(pval.matrix, order, plot=TRUE)
aggrPvals(pval.matrix, order, plot=TRUE)
pval.matrix |
Numeric matrix of p-values, columns represent different sets of p-values |
order |
Numeric constant, the order statistic that is used for the aggregation. |
plot |
Boolean value whether to plot p-value distributions. |
Aggregated p-value of the given order.
Daniela Beisser
data(pvaluesExample) aggrPvals(pval.matrix=pvaluesExample, order=2)
data(pvaluesExample) aggrPvals(pval.matrix=pvaluesExample, order=2)
The function fits a beta-uniform mixture model to a given p-value distribution.
bumOptim(x, starts=1, labels=NULL)
bumOptim(x, starts=1, labels=NULL)
x |
Numerical vector of p-values, has to be named with the gene names or the gene names can be given in the labels paramater. |
starts |
Number of start points for the optimization. |
labels |
Gene names for the p-values. |
List of class fb with the following elements:
lambda |
Fitted parameter lambda for the beta-uniform mixture model. |
a |
Fitted parameter a for the beta-uniform mixture model. |
negLL |
Negative log-likelihood. |
pvalues |
P-value vector. |
Marcus Dittrich and Daniela Beisser
M. T. Dittrich, G. W. Klau, A. Rosenwald, T. Dandekar, T. Mueller (2008) Identifying functional modules in protein-protein interaction networks: an integrated exact approach. (ISMB2008) Bioinformatics, 24: 13. i223-i231 Jul.
S. Pounds, S.W. Morris (2003) Estimating the occurrence of false positives and false negatives in microarray studies by approximating and partitioning the empirical distribution of p-values. Bioinformatics, 19(10): 1236-1242.
fitBumModel
, plot.bum
, hist.bum
data(pvaluesExample) pvals <- pvaluesExample[,1] bum <- bumOptim(x=pvals, starts=10) bum
data(pvaluesExample) pvals <- pvaluesExample[,1] bum <- bumOptim(x=pvals, starts=10) bum
The function compares the following parameters of two networks: diameter, average degree, degree exponent, average path length and plots the cumulative degree distributions. The networks have to be connected components.
compareNetworks(network1, network2, plot=TRUE)
compareNetworks(network1, network2, plot=TRUE)
network1 |
Network graphNEL or igraph format. |
network2 |
Second network in graphNEL or igraph format, or subnetwork drawn from first network. |
plot |
Boolean value, whether to plot the cumulative degree distributions. |
A vector of network parameters is returned:
diam.network1 |
Network diameter |
diam.network2 |
Diameter of the subnetwork |
av.degree.network1 |
Average degree of the network |
av.degree.network2 |
Average degree of the subnetwork |
degree.exponent.network1 |
Degree exponent of the network |
degree.exponent.network2 |
Degree exponent of the subnetwork |
av.path.length.network1 |
Average path lenght of the network |
av.path.length.network2 |
Average path length of the subnetwork |
Daniela Beisser
library(DLBCL) data(interactome) subnet1 <- largestComp(subNetwork(nodes(interactome)[1:100], interactome)) subnet2 <- largestComp(subNetwork(nodes(interactome)[101:200], interactome)) compareNetworks(network1=subnet1, network2=subnet2)
library(DLBCL) data(interactome) subnet1 <- largestComp(subNetwork(nodes(interactome)[1:100], interactome)) subnet2 <- largestComp(subNetwork(nodes(interactome)[101:200], interactome)) compareNetworks(network1=subnet1, network2=subnet2)
The function calculates consensus scores for a network, given a list of replicate modules.
consensusScores(modules, network, ro=length(modules)/2)
consensusScores(modules, network, ro=length(modules)/2)
modules |
Calculated modules from pseudo-replicates of expression values in igraph or graphNEL format. |
network |
Interaction network, which shoupld be scores. In igraph or graphNEL format |
ro |
Threshold which is subtracted from the scores to obtain positive and negative value. The default value is half of the number of replicates. |
A result list is returned, consisting of:
N.scores |
Numerical vector node scores. |
E.scores |
Numerical vector edge scores. |
N.frequencies |
Numerical vector node frequencies from the replicate modules. |
E.frequencies |
Numerical vector edge frequencies from the replicate modules. |
Daniela Beisser
library(DLBCL) data(interactome) network <- interactome # precomputed Heinz modules from pseudo-replicates ## Not run: lib <- file.path(.path.package("BioNet"), "extdata") modules <- readHeinzGraph(node.file=file.path(datadir, "ALL_n_resample.txt.0.hnz"), network=network) cons.scores <- consensusScores(modules, network) ## End(Not run)
library(DLBCL) data(interactome) network <- interactome # precomputed Heinz modules from pseudo-replicates ## Not run: lib <- file.path(.path.package("BioNet"), "extdata") modules <- readHeinzGraph(node.file=file.path(datadir, "ALL_n_resample.txt.0.hnz"), network=network) cons.scores <- consensusScores(modules, network) ## End(Not run)
Function to compute the density of the beta-uniform mixture model.
fbum(x, lambda, a)
fbum(x, lambda, a)
x |
A numeric value. |
lambda |
Parameter lambda, mixture parameter, proportion of uniform component |
a |
Parameter a, shape parameter of beta component |
Value of the density of the bum distribution for x.
Marcus Dittrich
S. Pounds, S.W. Morris (2003) Estimating the occurrence of false positives and false negatives in microarray studies by approximating and partitioning the empirical distribution of p-values. Bioinformatics, 19(10): 1236-1242.
y <- fbum(x=0.5, lambda=0.1, a=0.1) y
y <- fbum(x=0.5, lambda=0.1, a=0.1) y
The function calculates the log likelihood of the BUM model.
fbumLL(parms, x)
fbumLL(parms, x)
parms |
Vector of parameters; lambda and a. |
x |
Numerical vector of p-values. |
Log likelihood.
Marcus Dittrich
data(pvaluesExample) pvals <- pvaluesExample[,1] bum.mle <- fitBumModel(pvals, plot=FALSE) fbumLL(parms=c(bum.mle$lambda, bum.mle$a), x=pvals)
data(pvaluesExample) pvals <- pvaluesExample[,1] bum.mle <- fitBumModel(pvals, plot=FALSE) fbumLL(parms=c(bum.mle$lambda, bum.mle$a), x=pvals)
The function calculates the p-value threshold tau for a given false discovery rate. Tau is used for the scoring function.
fdrThreshold(fdr, fb)
fdrThreshold(fdr, fb)
fdr |
False discovery rate. |
fb |
Model from the beta-uniform mixture fitting. |
P-value threshold tau.
Marcus Dittrich
S. Pounds, S.W. Morris (2003) Estimating the occurrence of false positives and false negatives in microarray studies by approximating and partitioning the empirical distribution of p-values. Bioinformatics, 19(10): 1236-1242.
data(pvaluesExample) pvals <- pvaluesExample[,1] bum.mle <- fitBumModel(pvals, plot=FALSE) tau <- fdrThreshold(fdr=0.001, fb=bum.mle) tau
data(pvaluesExample) pvals <- pvaluesExample[,1] bum.mle <- fitBumModel(pvals, plot=FALSE) tau <- fdrThreshold(fdr=0.001, fb=bum.mle) tau
The function fits a beta-uniform mixture model to a given p-value distribution. The BUM method was introduced by Stan Pounds and Steve Morris to model the p-value distribution as a signal-noise decompostion. The signal component is assumed to be B(a,1)-distributed, whereas the noise component is uniform-distributed under the null hypothesis.
fitBumModel(x, plot = TRUE, starts=10)
fitBumModel(x, plot = TRUE, starts=10)
x |
Numeric vector of p-values. |
plot |
Boolean value, whether to plot a histogram and qqplot of the p-values with the fitted model. |
starts |
Numeric value giving the number of starts for the optimization. |
Maximum likelihood estimator object for the fitted bum model. List of class fb with the following elements:
lambda |
Fitted parameter lambda for the beta-uniform mixture model. |
a |
Fitted parameter a for the beta-uniform mixture model. |
negLL |
Negative log-likelihood. |
pvalues |
P-value vector. |
Daniela Beisser
S. Pounds, S.W. Morris (2003) Estimating the occurrence of false positives and false negatives in microarray studies by approximating and partitioning the empirical distribution of p-values. Bioinformatics, 19(10): 1236-1242.
data(pvaluesExample) pvals <- pvaluesExample[,1] bum.mle <- fitBumModel(pvals, plot=TRUE) bum.mle
data(pvaluesExample) pvals <- pvaluesExample[,1] bum.mle <- fitBumModel(pvals, plot=TRUE) bum.mle
The function partitions the scores into scores for each subgraph of the network.
getCompScores(network, score)
getCompScores(network, score)
network |
A network in graphNEL or igraph format. |
score |
Vector of scores. |
A data frame with the components of the network and the score for each PPI identifier.
Marcus Dittrich
library(DLBCL) data(interactome) data(dataLym) # create random subgraph with 100 nodes and their direct neighbors nodes <- nodes(interactome)[sample(length(nodes(interactome)), 100)] subnet <- subNetwork(nodeList=nodes, network=interactome, neighbors="first") score <- dataLym$score001 names(score) <- dataLym$label getCompScores(score=score, network=subnet)
library(DLBCL) data(interactome) data(dataLym) # create random subgraph with 100 nodes and their direct neighbors nodes <- nodes(interactome)[sample(length(nodes(interactome)), 100)] subnet <- subNetwork(nodeList=nodes, network=interactome, neighbors="first") score <- dataLym$score001 names(score) <- dataLym$label getCompScores(score=score, network=subnet)
A network in graphNEL or igraph format is converted to an edgelist.
getEdgeList(network)
getEdgeList(network)
network |
Network in graphNEL or igraph format. |
A matrix whose columns represent the connected edges.
Marcus Dittrich
library(DLBCL) data(interactome) getEdgeList(interactome)[1:10,]
library(DLBCL) data(interactome) getEdgeList(interactome)[1:10,]
The function plots a histogram of the p-values together with the fitted bum-model.
## S3 method for class 'bum' hist(x, breaks=50, main="Histogram of p-values", xlab="P-values", ylab="Density", ...)
## S3 method for class 'bum' hist(x, breaks=50, main="Histogram of p-values", xlab="P-values", ylab="Density", ...)
x |
Maximum likelihood estimator object of the beta-uniform mixture fit. |
breaks |
Breaks for the histogram. |
main |
An overall title for the plot. |
xlab |
A title for the x axis. |
ylab |
A title for the y axis. |
... |
Other graphic parameters for the plot. |
Daniela Beisser
fitBumModel
, hist.bum
, bumOptim
data(pvaluesExample) pvals <- pvaluesExample[,1] mle <- fitBumModel(pvals, plot=FALSE) hist(mle)
data(pvaluesExample) pvals <- pvaluesExample[,1] mle <- fitBumModel(pvals, plot=FALSE) hist(mle)
The function extracts the largest component of a network.
largestComp(network)
largestComp(network)
network |
A graph in graphNEL or igraph format. |
A new graph object that represents the largest component of the given network.
Marcus Dittrich
library(DLBCL) data(interactome) interactome largestComp(interactome)
library(DLBCL) data(interactome) interactome largestComp(interactome)
The function extracts the component of the network with the largest score. All nodes have to exceed the given level for the score.
largestScoreComp(network, score, level=0)
largestScoreComp(network, score, level=0)
network |
Network in graphNEL or igraph format. |
score |
Vector of scores for the network. |
level |
Cut-off level for the score for the component. |
Subgraph of the network with a score larger than the given level.
Marcus Dittrich
library(DLBCL) data(interactome) data(dataLym) network <- rmSelfLoops(interactome) score <- dataLym$score001 names(score) <- dataLym$label lComp <- largestScoreComp(network=network, score=score, level=1) ## Not run: plotModule(lComp)
library(DLBCL) data(interactome) data(dataLym) network <- rmSelfLoops(interactome) score <- dataLym$score001 names(score) <- dataLym$label lComp <- largestScoreComp(network=network, score=score, level=1) ## Not run: plotModule(lComp)
The function loads a network from a Cytoscape sif file. Edge attributes are provided in the ea.file or vector of ea.files. The node attributes are provided the same way. For other formats see read.graph in the igraph package.
loadNetwork.sif(sif.file, na.file, ea.file, format=c("graphNEL", "igraph"), directed=FALSE)
loadNetwork.sif(sif.file, na.file, ea.file, format=c("graphNEL", "igraph"), directed=FALSE)
sif.file |
Cytoscape sif file, containing the network. |
na.file |
File or vector of file with Cytoscape node attibutes. |
ea.file |
File or vector of file with Cytoscape edge attibutes. |
format |
Format of output graph, either graphNEL or igraph. |
directed |
Boolean value for directed or undirected graph. |
Graph with loaded attributes.
Daniela Beisser
## Not run: lib <- file.path(.path.package("BioNet"), "extdata") # load interaction file, node attribute file with a node weight of 2 for each node and the edge attribute file with a edge weight of 1 for each edge network <- loadNetwork.sif(sif.file=file.path(lib,"cytoscape.sif"), na.file=file.path(lib,"n.weight.NA"), ea.file=file.path(lib,"weight.EA"), format="graphNEL", directed=FALSE); network; nodeData(network); edgeData(network); ## End(Not run)
## Not run: lib <- file.path(.path.package("BioNet"), "extdata") # load interaction file, node attribute file with a node weight of 2 for each node and the edge attribute file with a edge weight of 1 for each edge network <- loadNetwork.sif(sif.file=file.path(lib,"cytoscape.sif"), na.file=file.path(lib,"n.weight.NA"), ea.file=file.path(lib,"weight.EA"), format="graphNEL", directed=FALSE); network; nodeData(network); edgeData(network); ## End(Not run)
The function loads a network from a tabular format.
loadNetwork.tab(file, header=TRUE, directed=FALSE, format=c("graphNEL", "igraph"))
loadNetwork.tab(file, header=TRUE, directed=FALSE, format=c("graphNEL", "igraph"))
file |
File with network to load. |
header |
Booelan value whether to include header or not. |
directed |
Booelan value whether the network is to be directed or not. |
format |
Output format of the network, either graphNEL or igraph |
Marcus Dittrich
Function to create a graph in graphNEL or igraph format from a source and a target vector.
makeNetwork(source, target, edgemode="undirected", format=c("graphNEL", "igraph"))
makeNetwork(source, target, edgemode="undirected", format=c("graphNEL", "igraph"))
source |
Vector of source nodes. |
target |
Vector of corresponding target nodes. |
edgemode |
For an "undirected" or "directed" network. |
format |
Graph format, eiter graphNEL or igraph. |
A graph object.
Marcus Dittrich
source <- c("a", "b", "c", "d") target <- c("b", "c", "a", "a") graph <- makeNetwork(source, target, edgemode="undirected")
source <- c("a", "b", "c", "d") target <- c("b", "c", "a", "a") graph <- makeNetwork(source, target, edgemode="undirected")
The function selects for each gene the probeset with the highest variance and gets the PPI ID for each gene. The PPI identifier is: gene symbol(Entrez ID). Affymetrix identifiers are mapped to the ENTREZ ID.
mapByVar(exprSet, network=NULL, attr="geneID", ignoreAFFX=TRUE)
mapByVar(exprSet, network=NULL, attr="geneID", ignoreAFFX=TRUE)
exprSet |
Affymetrix ExpressionSet. |
network |
Network that is used to map the Affymetrix identifiers. |
attr |
The attribute of the network that is used to map the Affymetrix IDs. The IDs are mapped to the unique Entrez gene IDs, which are by default stored in the "geneID" attribute of the network. |
ignoreAFFX |
Boolean value, whether to ignore or leave AFFX control genes. |
Expression matrix with one gene (PPI ID) per probeset.
Daniela Beisser
## Not run: library(ALL); data(ALL); mapped.e.set <- mapByVar(ALL); mapped.e.set[1:10,]; ## End(Not run)
## Not run: library(ALL); data(ALL); mapped.e.set <- mapByVar(ALL); mapped.e.set[1:10,]; ## End(Not run)
Function to permutate node labels of a given network.
permutateNodes(network)
permutateNodes(network)
network |
Network in graphNEL or igraph format. |
Network with permutated labels.
Marcus Dittrich
library(DLBCL) data(interactome) # remove self-loops before permutating the labels interactome <- rmSelfLoops(interactome) perm.net <- permutateNodes(interactome) perm.net
library(DLBCL) data(interactome) # remove self-loops before permutating the labels interactome <- rmSelfLoops(interactome) perm.net <- permutateNodes(interactome) perm.net
The function calculates the upper bound pi for the fraction of noise.
piUpper(fb)
piUpper(fb)
fb |
Fitted bum model, list with parameters a and lambda. |
Numerical value for the upper bound pi.
Marcus Dittrich
data(pvaluesExample) pvals <- pvaluesExample[,1] bum <- bumOptim(pvals, starts=10) piUpper(fb=bum)
data(pvaluesExample) pvals <- pvaluesExample[,1] bum <- bumOptim(pvals, starts=10) piUpper(fb=bum)
The function plot the theoretical quantiles of the fitted bum model against the quantiles of the observed p-value distribution.
## S3 method for class 'bum' plot(x, main="QQ-Plot", xlab="Estimated p-value", ylab="Observed p-value", ...)
## S3 method for class 'bum' plot(x, main="QQ-Plot", xlab="Estimated p-value", ylab="Observed p-value", ...)
x |
Maximum likelihood estimation object of the fitted bum model. |
main |
An overall title for the plot. |
xlab |
A title for the x axis. |
ylab |
A title for the y axis. |
... |
Other graphic parameters for the plot. |
Daniela Beisser
fitBumModel
, plot.bum
, bumOptim
data(pvaluesExample) pvals <- pvaluesExample[,1] mle <- fitBumModel(pvals, plot=FALSE) plot(mle)
data(pvaluesExample) pvals <- pvaluesExample[,1] mle <- fitBumModel(pvals, plot=FALSE) plot(mle)
The function plots a network from graphNEL or igraph format in 3D using a modified function from the package igraph and requires the package rgl which uses openGL. The 3D plot can be zoomed, rotated, shifted on the canvas. This function is just used to visualize the modules. For further plotting options use the rglplot function of the igraph package. If a score attribute is provided in the graph this will be used for the coloring of the nodes. Otherwise a vector of values can be given by the diff.or.score argument. The vector has to contain positive and negative values, either scores or values for differential expression (fold changes). Labels for the nodes can be provided by the labels argument, otherwise it will be automatically looked for a geneSymbol attribute of the nodes.
plot3dModule(network, labels=NULL, windowSize = c(100,100,1500,1000), diff.or.scores=NULL, red=c("negative", "positive"), ...)
plot3dModule(network, labels=NULL, windowSize = c(100,100,1500,1000), diff.or.scores=NULL, red=c("negative", "positive"), ...)
network |
Network in graphNEL or igraph format. |
labels |
Labels for the nodes of the network. Otherwise it will be automatically looked for a geneSymbol attribute of the nodes. |
windowSize |
Numerical vector of size four to set the size of the rgl device. |
diff.or.scores |
Named numerical vector of differential expression (fold changes) or scores of the nodes in the network. These will be used for node coloring. Otherwise a score attribute of the nodes will be automatically used. |
red |
Either "negative" or "positive", to specify which values are to be colored red in the plot. |
... |
Other graphic parameters for the plot. |
Daniela Beisser
library(DLBCL) data(interactome) data(dataLym) interactome <- subNetwork(dataLym$label, interactome) interactome <- rmSelfLoops(interactome) fchange <- dataLym$diff names(fchange) <- dataLym$label subnet <- largestComp(subNetwork(nodes(interactome)[1:100], interactome)) diff <- fchange[nodes(subnet)] ## Not run: library(rgl); plot3dModule(network=subnet, diff.or.score=diff) ## End(Not run)
library(DLBCL) data(interactome) data(dataLym) interactome <- subNetwork(dataLym$label, interactome) interactome <- rmSelfLoops(interactome) fchange <- dataLym$diff names(fchange) <- dataLym$label subnet <- largestComp(subNetwork(nodes(interactome)[1:100], interactome)) diff <- fchange[nodes(subnet)] ## Not run: library(rgl); plot3dModule(network=subnet, diff.or.score=diff) ## End(Not run)
The function plots the log likelihood surface for all a and lambda parameter of the beta-uniform mixture model.
plotLLSurface(x, opt=NULL, main="Log-Likelihood Surface", color.palette = heat.colors, nlevels = 32)
plotLLSurface(x, opt=NULL, main="Log-Likelihood Surface", color.palette = heat.colors, nlevels = 32)
x |
Numeric vector of p-values. |
opt |
List of optimal parameters for a and lambda from the beta-uniform mixture model. |
main |
The overall title of the plot. |
color.palette |
Color scheme of the image plot. |
nlevels |
Number of color levels. |
Marcus Dittrich
library(DLBCL) data(dataLym) pvals <- dataLym$t.pval names(pvals) <- dataLym$label mle <- fitBumModel(pvals, plot=FALSE) plotLLSurface(x=pvals, opt=mle)
library(DLBCL) data(dataLym) pvals <- dataLym$t.pval names(pvals) <- dataLym$label mle <- fitBumModel(pvals, plot=FALSE) plotLLSurface(x=pvals, opt=mle)
The function plots a network from graphNEL or igraph format, adapted from an igraph plotting function. It is just used to visualize the modules. For further plotting options use the plot.igraph function of the igraph package. The shapes of the nodes can be changed according to the scores argument, then negative scores appear squared. The color of the nodes can be changed according to the diff.expr argument. Negative values lead to green nodes, positive values are colored in red. If the vectors are not provided, it will be automatically looked for nodes attributes with the name score and diff.expr.
plotModule(network, layout=layout.fruchterman.reingold, labels=NULL, diff.expr=NULL, scores=NULL, main=NULL, vertex.size=NULL, ...)
plotModule(network, layout=layout.fruchterman.reingold, labels=NULL, diff.expr=NULL, scores=NULL, main=NULL, vertex.size=NULL, ...)
network |
Network in graphNEL or igraph format. |
layout |
Layout algorithm, e.g. layout.fruchterman.reingold or layout.kamada.kawai. |
labels |
Labels for the nodes of the network. |
diff.expr |
Named numerical vector of differential expression (fold changes) of the nodes in the network. These will be used for coloring of the nodes. It will be automatically looked for nodes attribute with the name diff.expr, if the argument is null. |
scores |
Named numerical vector of scores of the nodes in the network. These will be used for the shape of the nodes. It will be automatically looked for nodes attribute with the name score, if the argument is null. |
main |
Main title of the plot. |
vertex.size |
Numerical value or verctor for the size of the vertices. |
... |
Other graphic parameters for the plot. |
Marcus Dittrich and Daniela Beisser
library(DLBCL) data(dataLym) data(interactome) interactome <- subNetwork(dataLym$label, interactome) interactome <- rmSelfLoops(interactome) fchange <- dataLym$diff names(fchange) <- dataLym$label # create random subnetwork subnet <- largestComp(subNetwork(nodes(interactome)[1:100], interactome)) fchange <- fchange[nodes(subnet)] # color random subnetwork by the fold change ## Not run: plotModule(network=subnet, diff.expr=fchange)
library(DLBCL) data(dataLym) data(interactome) interactome <- subNetwork(dataLym$label, interactome) interactome <- rmSelfLoops(interactome) fchange <- dataLym$diff names(fchange) <- dataLym$label # create random subnetwork subnet <- largestComp(subNetwork(nodes(interactome)[1:100], interactome)) fchange <- fchange[nodes(subnet)] # color random subnetwork by the fold change ## Not run: plotModule(network=subnet, diff.expr=fchange)
The function prints information about the bum model.
## S3 method for class 'bum' print(x, ...)
## S3 method for class 'bum' print(x, ...)
x |
Maximum likelihood estimator object of the beta-uniform mixture fit. |
... |
Other graphic parameters for print. |
Marcus Dittrich
data(pvaluesExample) pvals <- pvaluesExample[,1] mle <- fitBumModel(pvals, plot=FALSE) print(mle)
data(pvaluesExample) pvals <- pvaluesExample[,1] mle <- fitBumModel(pvals, plot=FALSE) print(mle)
Data example consisting of a matrix of p-values. Each gene has two corresponding p-values.
These p-values can be aggregated into a p-value of p-values by the method aggrPvals
.
data(pvaluesExample)
data(pvaluesExample)
data(pvaluesExample) pvaluesExample[1:10,]
data(pvaluesExample) pvaluesExample[1:10,]
Function to convert the HEINZ output to a graph object, or if the output is in matrix form, it will create a list of graphs. The function needs the node and the original network, from which the module is calculated.
readHeinzGraph(node.file, network, format=c("graphNEL", "igraph"))
readHeinzGraph(node.file, network, format=c("graphNEL", "igraph"))
node.file |
Heinz node output file. |
network |
Original network from which Heinz input was created. |
format |
Graph format of output, either igraph or graphNEL. |
Graph object.
Daniela Beisser
library(DLBCL) data(interactome) # precomputed Heinz output files ## Not run: lib <- file.path(path.package("BioNet"), "extdata") module <- readHeinzGraph(node.file=file.path(lib, "lymphoma_nodes_001.txt.0.hnz"), network=interactome, format="graphNEL"); plotModule(module); ## End(Not run)
library(DLBCL) data(interactome) # precomputed Heinz output files ## Not run: lib <- file.path(path.package("BioNet"), "extdata") module <- readHeinzGraph(node.file=file.path(lib, "lymphoma_nodes_001.txt.0.hnz"), network=interactome, format="graphNEL"); plotModule(module); ## End(Not run)
Converts the HEINZ output to a tree in graph format. If the output is in matrix form, it will create a list of graphs. The function needs the node and edge file and the original network from which the module is calculated.
readHeinzTree(node.file, edge.file, network, format=c("graphNEL", "igraph"))
readHeinzTree(node.file, edge.file, network, format=c("graphNEL", "igraph"))
node.file |
Heinz node output file. |
edge.file |
Heinz edge output file. |
network |
Original network from which Heinz input was created. |
format |
Output format of the graph, either igraph or graphNEL. |
A graph object.
Daniela Beisser
library(DLBCL) data(interactome) # precomputed Heinz output files ## Not run: lib <- file.path(.path.package("BioNet"), "extdata") module <- readHeinzTree(node.file=file.path(lib, "lymphoma_nodes_001.txt.0.hnz"), edge.file=file.path(lib, "lymphoma_edges_001.txt.0.hnz"), network=interactome, format="graphNEL"); plotModule(module); ## End(Not run)
library(DLBCL) data(interactome) # precomputed Heinz output files ## Not run: lib <- file.path(.path.package("BioNet"), "extdata") module <- readHeinzTree(node.file=file.path(lib, "lymphoma_nodes_001.txt.0.hnz"), edge.file=file.path(lib, "lymphoma_edges_001.txt.0.hnz"), network=interactome, format="graphNEL"); plotModule(module); ## End(Not run)
The function uses a 50% jackknife resampling to calculate a pseudo-replicate of the expression matrix. The resampled expression values are used thereupon to calculate p-values for the differential expression between the given groups. Only two-group comparisons are allowed for the perfomed t-test.
resamplingPvalues(exprMat, groups, alternative = c("two.sided", "less", "greater"), resampleMat=FALSE)
resamplingPvalues(exprMat, groups, alternative = c("two.sided", "less", "greater"), resampleMat=FALSE)
exprMat |
Matrix with microarray expression values. |
groups |
Factors for two groups that are tested for differential expression. |
alternative |
Testing alternatives for the t-test: "two.sided", "less" or "greater". |
resampleMat |
Boolean value, whether to retrieve the matrix of jacknife resamples or not. |
A result list is returned, consisting of:
p.values |
VNumerical vector of p-values. |
resampleMat |
Matrix of resampled expression values. |
Daniela Beisser
library(ALL) data(ALL) mat <- exprs(ALL) groups <- factor(c(rep("A", 64), rep("B", 64))) results <- resamplingPvalues(mat, groups, alternative="greater")
library(ALL) data(ALL) mat <- exprs(ALL) groups <- factor(c(rep("A", 64), rep("B", 64))) results <- resamplingPvalues(mat, groups, alternative="greater")
The function removes self-loops, edges that start and end in the same node, from the network.
rmSelfLoops(network)
rmSelfLoops(network)
network |
A graph object, either in graphNEL or igraph format. |
The graph with the removed edges.
Marcus Dittrich
graph <- makeNetwork(c("a","b","c","d","e","a"), c("b","c","d","e","e","e")) graph2 <- rmSelfLoops(graph) edges(graph) edges(graph2)
graph <- makeNetwork(c("a","b","c","d","e","a"), c("b","c","d","e","e","e")) graph2 <- rmSelfLoops(graph) edges(graph) edges(graph2)
The function uses an heuristic approach to calculate the maximum scoring subnetwork. Based on the given network and scores the positive nodes are in the first step aggregated to meta-nodes between which minimum spanning trees are calculated. In regard to this, shortest paths yield the approximated maximum scoring subnetwork. This function can be used if a CPLEX license is not available to calculate the optimal solution.
runFastHeinz(network, scores)
runFastHeinz(network, scores)
network |
A graph in igraph or graphNEL format. |
scores |
A named vector, containing the scores for the nodes of the network. All nodes need to be scored in order to run the algorithm. |
A subnetwork in the input network format.
Daniela Beisser
writeHeinzEdges
, writeHeinzNodes
, readHeinzTree
, readHeinzGraph
, runHeinz
library(DLBCL) # load p-values data(dataLym) # load graph data(interactome) # get induced subnetwork for all genes contained on the chip interactome <- subNetwork(dataLym$label, interactome) p.values <- dataLym$t.pval names(p.values) <- dataLym$label bum <- fitBumModel(p.values, plot=TRUE) scores <- scoreNodes(network=interactome, fb=bum, fdr=0.0001) module <- runFastHeinz(network=interactome, scores=scores) ## Not run: plotModule(module)
library(DLBCL) # load p-values data(dataLym) # load graph data(interactome) # get induced subnetwork for all genes contained on the chip interactome <- subNetwork(dataLym$label, interactome) p.values <- dataLym$t.pval names(p.values) <- dataLym$label bum <- fitBumModel(p.values, plot=TRUE) scores <- scoreNodes(network=interactome, fb=bum, fdr=0.0001) module <- runFastHeinz(network=interactome, scores=scores) ## Not run: plotModule(module)
The function starts HEINZ from command line. The HEINZ folder has to include the heinz.py python script and the dhea file. CPLEX has to be installed and accessible from the computer R runs on.
runHeinz(heinz.folder="", heinz.e.file, heinz.n.file, N=TRUE, E=FALSE, diff=-1, n=1)
runHeinz(heinz.folder="", heinz.e.file, heinz.n.file, N=TRUE, E=FALSE, diff=-1, n=1)
heinz.folder |
The folder which contains the heinz.py python script and the dhea file. |
heinz.e.file |
The HEINZ edge input file. See |
heinz.n.file |
The HEINZ node input file. See |
N |
Boolean value, whether to run HEINZ on nodes. |
E |
Boolean value, whether to run HEINZ on edges. HEINZ can run on both with N and E set to TRUE. |
diff |
Difference of suboptimal solutions to optimal solution in haming distance in percent. Parameter is set to -1 for optimal solution. |
n |
Number of optimal and suboptimal solutions, the standard n=1 delivers only the optimal solution. |
This function starts the integer linear programming algorithm to calculate the optimal scoring subnetwork.
The algorithm might be started in the command line when the CPLEX is installed on another machine. To start
it from command line use: heinz.py -e edge.file.txt -n node.file.txt -E False/True -N False/True.
The results can be loaded with readHeinzTree
, readHeinzGraph
as a graph object.
Daniela Beisser
M. T. Dittrich, G. W. Klau, A. Rosenwald, T. Dandekar, T. Mueller (2008) Identifying functional modules in protein-protein interaction networks: an integrated exact approach. (ISMB2008) Bioinformatics, 24: 13. i223-i231 Jul.
writeHeinzEdges
, writeHeinzNodes
, readHeinzTree
, readHeinzGraph
The function saves a 3D plot of a network to file, therefore it requires the plot to be open. A screenshot of the 3D plot can be saved in "pdf" format. Background of the device is changed to white for plotting. The screenshot can take several seconds for large plots.
save3dModule(file)
save3dModule(file)
file |
File to save to. |
Daniela Beisser
library(DLBCL) data(dataLym) data(interactome) interactome <- subNetwork(dataLym$label, interactome) fchange <- dataLym$diff names(fchange) <- dataLym$label subnet <- largestComp(subNetwork(nodes(interactome)[1:100], interactome)) diff <- fchange[nodes(subnet)] ## Not run: library(rgl); plot3dModule(network=subnet, diff.or.score=diff); save3dModule(file="test") ## End(Not run)
library(DLBCL) data(dataLym) data(interactome) interactome <- subNetwork(dataLym$label, interactome) fchange <- dataLym$diff names(fchange) <- dataLym$label subnet <- largestComp(subNetwork(nodes(interactome)[1:100], interactome)) diff <- fchange[nodes(subnet)] ## Not run: library(rgl); plot3dModule(network=subnet, diff.or.score=diff); save3dModule(file="test") ## End(Not run)
The function saves a graph in a Cytoscape readable format: either in XGMML format, or as two tables, one for the nodes with attributes and one for the edges with attributes, or as .sif file. Or other standard formats like tab separated, .tgf, .net
saveNetwork(network, name="network", file, type=c("table", "XGMML", "sif", "tab", "tgf", "net"))
saveNetwork(network, name="network", file, type=c("table", "XGMML", "sif", "tab", "tgf", "net"))
network |
Network to save. |
name |
Name of the network, only needed for the XGMML format. |
file |
File to save to. |
type |
Type in which graph shall be saved. |
The format types are "XGMML", "table", "sif", "tab", "tgf" and "net". XGMML (eXtensible Graph Markup and Modeling Language) is an XML format based on GML which is used for graph description. Edges, nodes and their affiliated attributes are all saved in one file. In the table format two tables are created, one for the nodes with attributes and one for the edges with attributes. The sif format creates a .sif file for the network and an node attribute (.NA) or edge attribute (.EA) for each attribute. The name of the attribute is the filename. Tab writes only the edges of the network in a tabular format. Tgf save the network to simple .tgf format. The net format writes a Pajek readable file of the network and the ET type saves the edge tags to file.
Daniela Beisser and Marcus Dittrich
library(DLBCL) #create small network library(igraph) data(interactome) interactome <- igraph.from.graphNEL(interactome) small.net <- subNetwork(V(interactome)[1:16]$name, interactome) E(small.net)$e.weight <- rep(1,length(E(small.net))) V(small.net)$n.weight <- rep(2,length(V(small.net))) summary(small.net) ## Not run: saveNetwork(small.net, file="example_network", name="small.net", type="XGMML")
library(DLBCL) #create small network library(igraph) data(interactome) interactome <- igraph.from.graphNEL(interactome) small.net <- subNetwork(V(interactome)[1:16]$name, interactome) E(small.net)$e.weight <- rep(1,length(E(small.net))) V(small.net)$n.weight <- rep(2,length(V(small.net))) summary(small.net) ## Not run: saveNetwork(small.net, file="example_network", name="small.net", type="XGMML")
The function generates a dataframe for a given range of FDRs.
scanFDR(fb, fdr, labels=names(fb$pvalues))
scanFDR(fb, fdr, labels=names(fb$pvalues))
fb |
Fitted bum model. |
fdr |
Vector of FDRs. |
labels |
Data frame labels. |
Dataframe of scores for given p-values and a range of FDRs.
Marcus Dittrich
data(pvaluesExample) pvals <- pvaluesExample[,1] bum <- bumOptim(pvals, starts=10) scores <- scanFDR(fb=bum, fdr=c(0.1, 0.001, 0.0001)) scores[1:10,]
data(pvaluesExample) pvals <- pvaluesExample[,1] bum <- bumOptim(pvals, starts=10) scores <- scanFDR(fb=bum, fdr=c(0.1, 0.001, 0.0001)) scores[1:10,]
The function calculates a score for each gene with a given FDR from the fitted beta-uniform mixture model.
scoreFunction(fb, fdr=0.01)
scoreFunction(fb, fdr=0.01)
fb |
Model from the beta-uniform mixture fitting. |
fdr |
Numeric constant, from the false discovery rate a p-value threshold is calculated. P-values below this threshold are considered to be significant and will score positively, p-values a bove the threshold are supposed to arise from the null model. The FDR can be used to control the size of the maximum scoring subnetwork, by zooming in and out in the same region. |
Score vector for the given p-values.
Marcus Dittrich and Daniela Beisser
For details on the score calculation see: M. T. Dittrich, G. W. Klau, A. Rosenwald, T. Dandekar, T. Mueller (2008) Identifying functional modules in protein-protein interaction networks: an integrated exact approach. (ISMB2008) Bioinformatics, 24: 13. i223-i231 Jul.
data(pvaluesExample) pvals <- pvaluesExample[,1] bum.mle <- fitBumModel(pvals, plot=FALSE) scores <- scoreFunction(fdr=0.1, fb=bum.mle) scores
data(pvaluesExample) pvals <- pvaluesExample[,1] bum.mle <- fitBumModel(pvals, plot=FALSE) scores <- scoreFunction(fdr=0.1, fb=bum.mle) scores
The function derives scores from the p-values of the nodes of a network.
scoreNodes(network, fb, fdr=0.05)
scoreNodes(network, fb, fdr=0.05)
network |
A network in graphNEL or igraph format. |
fb |
Fitted bum model. |
fdr |
False discovery rate. |
Ordered score vector for the nodes of the network.
Marcus Dittrich
library(DLBCL) # load p-values data(dataLym) # load graph data(interactome) # get induced subnetwork for all genes contained on the chip chipGraph <- subNetwork(dataLym$label, interactome) p.values <- dataLym$t.pval names(p.values) <- dataLym$label bum <- fitBumModel(p.values, plot=TRUE) scoreNodes(network=chipGraph, fb=bum, fdr=0.001)
library(DLBCL) # load p-values data(dataLym) # load graph data(interactome) # get induced subnetwork for all genes contained on the chip chipGraph <- subNetwork(dataLym$label, interactome) p.values <- dataLym$t.pval names(p.values) <- dataLym$label bum <- fitBumModel(p.values, plot=TRUE) scoreNodes(network=chipGraph, fb=bum, fdr=0.001)
Function to change score offset from FDR1 to FDR2.
scoreOffset(fb, fdr1, fdr2)
scoreOffset(fb, fdr1, fdr2)
fb |
Model from the beta-uniform mixture fitting. |
fdr1 |
First false discovery rate. |
fdr2 |
Second false discovery rate. |
Offset for the score of the second FDR.
Marcus Dittrich
data(pvaluesExample) pvals <- pvaluesExample[,1] bum <- bumOptim(pvals, starts=10) scoreOffset(bum, fdr1=0.001, fdr2=0.000001)
data(pvaluesExample) pvals <- pvaluesExample[,1] bum <- bumOptim(pvals, starts=10) scoreOffset(bum, fdr1=0.001, fdr2=0.000001)
Function to get a sorted edgelist where the source protein is alphabetically smaller than the target protein from an undirected network.
sortedEdgeList(network)
sortedEdgeList(network)
network |
Undirected network in igraph or graphNEL format. |
Vector of sorted edges, where the source protein is alphabetically smaller than the target protein.
Daniela Beisser
library(DLBCL) data(interactome) E.list <- sortedEdgeList(interactome)
library(DLBCL) data(interactome) E.list <- sortedEdgeList(interactome)
The function creates a subgraph with the nodes given in the nodeList or for these nodes including their direct neighbors.
subNetwork(nodeList, network, neighbors=c("none", "first"))
subNetwork(nodeList, network, neighbors=c("none", "first"))
nodeList |
Character vector of nodes, contained in the subgraph. |
network |
Graph that is used for subgraph extraction. |
neighbors |
Neighborhood, that is chosen for the subgraph extraction. "none" are only the selected nodes, "first" includes the direct neighbors of the selected nodes. |
A graph object.
Marcus Dittrich
library(igraph) el <- cbind(c("a", "b", "c", "d", "e", "f", "d"), c("b", "c", "d", "e", "f", "a", "b")) graph <- graph.edgelist(el, directed=TRUE) node.list <- c("a", "b", "c") graph2 <- subNetwork(nodeList=node.list, network=graph) ## Not run: par(mfrow=c(1,2)); plotModule(graph); plotModule(graph2) ## End(Not run) # or in graphNEL format: graph3 <- igraph.to.graphNEL(graph) graph4 <- subNetwork(nodeList=node.list, network=graph3) graph3 graph4
library(igraph) el <- cbind(c("a", "b", "c", "d", "e", "f", "d"), c("b", "c", "d", "e", "f", "a", "b")) graph <- graph.edgelist(el, directed=TRUE) node.list <- c("a", "b", "c") graph2 <- subNetwork(nodeList=node.list, network=graph) ## Not run: par(mfrow=c(1,2)); plotModule(graph); plotModule(graph2) ## End(Not run) # or in graphNEL format: graph3 <- igraph.to.graphNEL(graph) graph4 <- subNetwork(nodeList=node.list, network=graph3) graph3 graph4
The function summarizes information about the bum model.
## S3 method for class 'bum' summary(object, ...)
## S3 method for class 'bum' summary(object, ...)
object |
Maximum likelihood estimator object of the beta-uniform mixture fit. |
... |
Other graphic parameters for summary. |
Daniela Beisser
data(pvaluesExample) pvals <- pvaluesExample[,1] mle <- fitBumModel(pvals, plot=FALSE) summary(mle)
data(pvaluesExample) pvals <- pvaluesExample[,1] mle <- fitBumModel(pvals, plot=FALSE) summary(mle)
Function to write the input files with the node and edge scores for HEINZ. These files are used
to calculate the maximum scoring subnetwork of the graph. The node scores are matched by their names to the nodes of the network, therefore if
nodes.scores are provided as a vector or matrix, the vector has to be named, respectively the matrix has to be provided
with rownames. If the network contains more nodes than the score vector, the nodes without a score are scored
with the average over all nodes. If the nodes should not be scored and used for the calculation of the
maximum scoring subnetwork, draw a subnetwork (subNetwork
) first and use this for the argument network.
The edge scores can be provided as a vector or matrix as the edge.scores argument. If no scores are provided in the arguments, but
the use.node.scores or use.edge.scores argument is set to TRUE, it will be automatically looked for the "score" attribute of the nodes
and edges of the network.
writeHeinz(network, file, node.scores=0, edge.scores=0, use.node.score=FALSE, use.edge.score=FALSE)
writeHeinz(network, file, node.scores=0, edge.scores=0, use.node.score=FALSE, use.edge.score=FALSE)
network |
Network from which to calculate the maximum scoring subnetwork. |
file |
File to write to. |
node.scores |
Numeric vector or matrix of scores for the nodes of the network. Names of the vector or rownames of the matrix have to correspond to the PPI identifiers of the network. The scores can also be used from the node attribute "score", given one score for each node. |
edge.scores |
Numeric vector of scores for the edges of the network. Edge scores have to be given in the order of the edges in the network. It is better to append the edge scores as the edge attribute "score" to the network: V(network)\$score and set use.scores to TRUE. |
use.node.score |
Boolean value, whether to use the node attribute "score" in the network as node scores. |
use.edge.score |
Boolean value, whether to use the edge attribute "score" in the network as edge scores. |
Daniela Beisser
writeHeinzNodes
and writeHeinzEdges
library(DLBCL) # use Lymphoma data and graph to find module data(interactome) data(dataLym) # get induced subnetwork for all genes contained on the chip chipGraph <- subNetwork(dataLym$label, interactome) score <- dataLym$score001 names(score) <- dataLym$label ## Not run: writeHeinz(network=chipGraph, file="lymphoma_001", node.scores=score, edge.scores=0)
library(DLBCL) # use Lymphoma data and graph to find module data(interactome) data(dataLym) # get induced subnetwork for all genes contained on the chip chipGraph <- subNetwork(dataLym$label, interactome) score <- dataLym$score001 names(score) <- dataLym$label ## Not run: writeHeinz(network=chipGraph, file="lymphoma_001", node.scores=score, edge.scores=0)
Function to write an input file for HEINZ with edge scores. If no edge scores are used, they are set to 0. In order to run HEINZ, a node input and edge input file are needed.
writeHeinzEdges(network, file, edge.scores=0, use.score=FALSE)
writeHeinzEdges(network, file, edge.scores=0, use.score=FALSE)
network |
Network from which to calculate the maximum scoring subnetwork. |
file |
File to write to. |
edge.scores |
Numeric vector of scores for the edges of the network. Edge scores have to be given in the order of the edges in the network. It is better to append the edge scores as the edge attribute "score" to the network: V(network)\$score and set use.score to TRUE. |
use.score |
Boolean value, whether to use the edge attribute "score" in the network as edge scores. |
Daniela Beisser
writeHeinzNodes
and writeHeinz
library(DLBCL) # use Lymphoma data and graph to find module data(interactome) data(dataLym) # get induced subnetwork for all genes contained on the chip chipGraph <- subNetwork(dataLym$label, interactome) # remove self loops graph <- rmSelfLoops(chipGraph) ## Not run: writeHeinzEdges(network=graph, file="lymphoma_edges_001", use.score=FALSE) score <- dataLym$score001 names(score) <- dataLym$label ## Not run: writeHeinzNodes(network=graph, file="lymphoma_nodes_001", node.scores=score) # write another edge file with edge scores library(igraph) data(interactome) interactome <- igraph.from.graphNEL(interactome) small.net <- subNetwork(V(interactome)[1:16]$name, interactome) scores <- c(1:length(E(small.net))) E(small.net)$score <- scores ## Not run: writeHeinzEdges(network=small.net, file="test_edges", use.score=TRUE)
library(DLBCL) # use Lymphoma data and graph to find module data(interactome) data(dataLym) # get induced subnetwork for all genes contained on the chip chipGraph <- subNetwork(dataLym$label, interactome) # remove self loops graph <- rmSelfLoops(chipGraph) ## Not run: writeHeinzEdges(network=graph, file="lymphoma_edges_001", use.score=FALSE) score <- dataLym$score001 names(score) <- dataLym$label ## Not run: writeHeinzNodes(network=graph, file="lymphoma_nodes_001", node.scores=score) # write another edge file with edge scores library(igraph) data(interactome) interactome <- igraph.from.graphNEL(interactome) small.net <- subNetwork(V(interactome)[1:16]$name, interactome) scores <- c(1:length(E(small.net))) E(small.net)$score <- scores ## Not run: writeHeinzEdges(network=small.net, file="test_edges", use.score=TRUE)
Function to write an input file with the node scores for HEINZ. This file is used
together with the edge input file to calculate the maximum scoring subnetwork of the
graph. The scores are matched by their names to the nodes of the network, therefore if
nodes.scores are provided as a vector or matrix, the vector has to be named, respectively the matrix has to be provided
with rownames.
If the network contains more nodes than the score vector, the nodes without a score are scored
with the average over all nodes. If the nodes should not be scored and used for the calculation of the
maximum scoring subnetwork, draw a subnetwork subNetwork
first and use this for the argument network.
writeHeinzNodes(network, file, node.scores=0, use.score=FALSE)
writeHeinzNodes(network, file, node.scores=0, use.score=FALSE)
network |
Network from which to calculate the maximum scoring subnetwork. |
file |
File to write to. |
node.scores |
Numeric vector or matrix of scores for the nodes of the network. Names of the vector or rownames of the matrix have to correspond to the PPI identifiers of the network. The scores can also be used from the node attribute "score", given one score for each node. |
use.score |
Boolean value, whether to use the node attribute "score" in the network as node scores. |
Use scoreNodes
or scoreFunction
to derive scores from a vector of p-values.
Daniela Beisser
writeHeinzEdges
and writeHeinz
#create small network library(DLBCL) data(interactome) small.net <- subNetwork(nodes(interactome)[0:15], interactome) scores <- c(1:length(nodes(small.net))) names(scores) <- nodes(small.net) ## Not run: writeHeinzNodes(network=small.net, file="test_nodes", node.scores=scores) # use Lymphoma data and graph to find module library(DLBCL) data(interactome) data(dataLym) # get induced subnetwork for all genes contained on the chip chipGraph <- subNetwork(dataLym$label, interactome) ## Not run: writeHeinzEdges(network=chipGraph, file="lymphoma_edges_001", use.score=FALSE) score <- dataLym$score001 names(score) <- dataLym$label ## Not run: writeHeinzNodes(network=chipGraph, file="lymphoma_nodes_001", node.scores=score)
#create small network library(DLBCL) data(interactome) small.net <- subNetwork(nodes(interactome)[0:15], interactome) scores <- c(1:length(nodes(small.net))) names(scores) <- nodes(small.net) ## Not run: writeHeinzNodes(network=small.net, file="test_nodes", node.scores=scores) # use Lymphoma data and graph to find module library(DLBCL) data(interactome) data(dataLym) # get induced subnetwork for all genes contained on the chip chipGraph <- subNetwork(dataLym$label, interactome) ## Not run: writeHeinzEdges(network=chipGraph, file="lymphoma_edges_001", use.score=FALSE) score <- dataLym$score001 names(score) <- dataLym$label ## Not run: writeHeinzNodes(network=chipGraph, file="lymphoma_nodes_001", node.scores=score)