Title: | plot bayesian network inferred from gene expression data based on enrichment analysis results |
---|---|
Description: | This package provides the visualization of bayesian network inferred from gene expression data. The networks are based on enrichment analysis results inferred from packages including clusterProfiler and ReactomePA. The networks between pathways and genes inside the pathways can be inferred and visualized. |
Authors: | Noriaki Sato [cre, aut] |
Maintainer: | Noriaki Sato <[email protected]> |
License: | Artistic-2.0 |
Version: | 1.7.0 |
Built: | 2024-11-19 03:30:22 UTC |
Source: | https://github.com/bioc/CBNplot |
Plot gene relationship within the specified pathway
bngeneplot( results, exp, expSample = NULL, algo = "hc", R = 20, returnNet = FALSE, algorithm.args = NULL, bypassConverting = FALSE, edgeLink = FALSE, pathNum = NULL, convertSymbol = TRUE, expRow = "ENSEMBL", interactive = FALSE, cexCategory = 1, cl = NULL, showDir = FALSE, chooseDir = FALSE, scoreType = "bic-g", labelSize = 4, layout = "nicely", clusterAlpha = 0.2, strType = "normal", delZeroDegree = TRUE, otherVar = NULL, otherVarName = NULL, onlyDf = FALSE, disc = FALSE, tr = NULL, remainCont = NULL, sp = "hsapiens", compareRef = FALSE, compareRefType = "intersection", pathDb = "reactome", dep = NULL, depMeta = NULL, sizeDep = FALSE, showDepHist = TRUE, cellLineName = "5637_URINARY_TRACT", showLineage = FALSE, orgDb = org.Hs.eg.db, shadowText = TRUE, bgColor = "white", textColor = "black", strengthPlot = FALSE, nStrength = 10, strThresh = NULL, hub = NULL, seed = 1, useSiGN = FALSE )
bngeneplot( results, exp, expSample = NULL, algo = "hc", R = 20, returnNet = FALSE, algorithm.args = NULL, bypassConverting = FALSE, edgeLink = FALSE, pathNum = NULL, convertSymbol = TRUE, expRow = "ENSEMBL", interactive = FALSE, cexCategory = 1, cl = NULL, showDir = FALSE, chooseDir = FALSE, scoreType = "bic-g", labelSize = 4, layout = "nicely", clusterAlpha = 0.2, strType = "normal", delZeroDegree = TRUE, otherVar = NULL, otherVarName = NULL, onlyDf = FALSE, disc = FALSE, tr = NULL, remainCont = NULL, sp = "hsapiens", compareRef = FALSE, compareRefType = "intersection", pathDb = "reactome", dep = NULL, depMeta = NULL, sizeDep = FALSE, showDepHist = TRUE, cellLineName = "5637_URINARY_TRACT", showLineage = FALSE, orgDb = org.Hs.eg.db, shadowText = TRUE, bgColor = "white", textColor = "black", strengthPlot = FALSE, nStrength = 10, strThresh = NULL, hub = NULL, seed = 1, useSiGN = FALSE )
results |
the enrichment analysis result |
exp |
gene expression matrix |
expSample |
candidate samples to be included in the inference default to all |
algo |
structure learning method used in boot.strength() default to "hc" |
R |
the number of bootstrap |
returnNet |
whether to return the network |
algorithm.args |
parameters to pass to bnlearn structure learnng function |
bypassConverting |
bypass the symbol converting If you use custom annotation databases that does not have SYMBOL listed in keys. ID of rownames and those listed in EA result must be same. |
edgeLink |
use geom_edge_link() instead of geom_edge_diagonal() |
pathNum |
the pathway number (the number of row of the original result, ordered by p-value) |
convertSymbol |
whether the label of resulting network is converted to symbol, default to TRUE |
expRow |
the type of the identifier of rows of expression matrix |
interactive |
whether to use bnviewer (default to FALSE) |
cexCategory |
scaling factor of size of nodes |
cl |
cluster object from parallel::makeCluster() |
showDir |
show the confidence of direction of edges |
chooseDir |
if undirected edges are present, choose direction of edges (default: FALSE) |
scoreType |
score type to use on choosing direction |
labelSize |
the size of label of the nodes |
layout |
ggraph layout, default to "nicely" |
clusterAlpha |
if specified multiple pathways, the parameter is passed to geom_mark_hull() |
strType |
"normal" or "ms" for multiscale implementation |
delZeroDegree |
delete zero degree nodes |
otherVar |
other variables to be included in the inference |
otherVarName |
the names of other variables |
onlyDf |
return only data.frame used for inference |
disc |
discretize the expressoin data |
tr |
Specify data.frame if one needs to discretize as the same parametersas the other dataset |
remainCont |
Specify characters when perform discretization, if some columns are to be remain continuous |
sp |
query to graphite::pathways(), default to "hsapiens" |
compareRef |
whether compare to the reference network |
compareRefType |
"intersection" or "difference" |
pathDb |
query to graphite::pathways(), default to "reactome" |
dep |
the tibble storing dependency score from library depmap |
depMeta |
depmap::depmap_metadata(), needed for showLineage |
sizeDep |
whether to reflect DepMap score to the node size |
showDepHist |
whether to show depmap histogram |
cellLineName |
the cell line name to be included |
showLineage |
show the dependency score across the lineage |
orgDb |
perform clusterProfiler::setReadable based on this organism database |
shadowText |
whether to use shadow text for the better readability default: TRUE |
bgColor |
color for text background when shadowText is TRUE |
textColor |
color for text when shadowText is TRUE |
strengthPlot |
append the barplot depicting edges with high strength |
nStrength |
specify how many edges are included in the strength plot |
strThresh |
the threshold for strength |
hub |
visualize the genes with top-n hub scores |
seed |
A random seed to make the analysis reproducible, default is 1. |
useSiGN |
default to FALSE. For using SiGN-BN in the function in Windows 10/11, 1. Download the SiGN-BN HC+BS binary in WSL (https://sign.hgc.jp/signbn/download.html) 2. Set PATH to executable (sign.1.8.3) |
ggplot2 object
data("exampleEaRes");data("exampleGeneExp") res <- bngeneplot(results = exampleEaRes, exp = exampleGeneExp, pathNum = 1, R = 10, convertSymbol = TRUE, expRow = "ENSEMBL")
data("exampleEaRes");data("exampleGeneExp") res <- bngeneplot(results = exampleEaRes, exp = exampleGeneExp, pathNum = 1, R = 10, convertSymbol = TRUE, expRow = "ENSEMBL")
Plot gene relationship within the specified pathway using customized theme
bngeneplotCustom( results, exp, expSample = NULL, algo = "hc", R = 20, pathNum = NULL, convertSymbol = TRUE, expRow = "ENSEMBL", interactive = FALSE, cexCategory = 1, cl = NULL, showDir = FALSE, chooseDir = FALSE, algorithm.args = NULL, labelSize = 4, layout = "nicely", strType = "normal", returnNet = FALSE, otherVar = NULL, otherVarName = NULL, onlyDf = FALSE, disc = FALSE, tr = NULL, remainCont = NULL, dep = NULL, sizeDep = FALSE, orgDb = org.Hs.eg.db, bypassConverting = FALSE, edgeLink = FALSE, cellLineName = "5637_URINARY_TRACT", fontFamily = "sans", strengthPlot = FALSE, nStrength = 10, strThresh = NULL, hub = NULL, glowEdgeNum = NULL, nodePal = c("blue", "red"), edgePal = c("blue", "red"), textCol = "black", titleCol = "black", backCol = "white", barTextCol = "black", barPal = c("red", "blue"), barBackCol = "white", scoreType = "bic-g", barLegendKeyCol = "white", barAxisCol = "black", bg.colour = NULL, bg.r = 0.1, barPanelGridCol = "black", titleSize = 24, seed = 1 )
bngeneplotCustom( results, exp, expSample = NULL, algo = "hc", R = 20, pathNum = NULL, convertSymbol = TRUE, expRow = "ENSEMBL", interactive = FALSE, cexCategory = 1, cl = NULL, showDir = FALSE, chooseDir = FALSE, algorithm.args = NULL, labelSize = 4, layout = "nicely", strType = "normal", returnNet = FALSE, otherVar = NULL, otherVarName = NULL, onlyDf = FALSE, disc = FALSE, tr = NULL, remainCont = NULL, dep = NULL, sizeDep = FALSE, orgDb = org.Hs.eg.db, bypassConverting = FALSE, edgeLink = FALSE, cellLineName = "5637_URINARY_TRACT", fontFamily = "sans", strengthPlot = FALSE, nStrength = 10, strThresh = NULL, hub = NULL, glowEdgeNum = NULL, nodePal = c("blue", "red"), edgePal = c("blue", "red"), textCol = "black", titleCol = "black", backCol = "white", barTextCol = "black", barPal = c("red", "blue"), barBackCol = "white", scoreType = "bic-g", barLegendKeyCol = "white", barAxisCol = "black", bg.colour = NULL, bg.r = 0.1, barPanelGridCol = "black", titleSize = 24, seed = 1 )
results |
the enrichment analysis result |
exp |
gene expression matrix |
expSample |
candidate rows to be included in the inference default to all |
algo |
structure learning method used in boot.strength() default to "hc" |
R |
the number of bootstrap |
pathNum |
the pathway number (the number of row of the original result, ordered by p-value) |
convertSymbol |
whether the label of resulting network is converted to symbol, default to TRUE |
expRow |
the type of the identifier of rows of expression matrix |
interactive |
whether to use bnviewer (default to FALSE) |
cexCategory |
scaling factor of size of nodes |
cl |
cluster object from parallel::makeCluster() |
showDir |
show the confidence of direction of edges |
chooseDir |
if undirected edges are present, choose direction of edges |
algorithm.args |
parameters to pass to bnlearn structure learnng function |
labelSize |
the size of label of the nodes |
layout |
ggraph layout, default to "nicely" |
strType |
"normal" or "ms" for multiscale implementation |
returnNet |
whether to return the network |
otherVar |
other variables to be included in the inference |
otherVarName |
the names of other variables |
onlyDf |
return only data.frame used for inference |
disc |
discretize the expressoin data |
tr |
Specify data.frame if one needs to discretize as the same parameters as the other dataset |
remainCont |
Specify characters when perform discretization, if some columns are to be remain continuous |
dep |
the tibble storing dependency score from library depmap |
sizeDep |
whether to reflect DepMap score to the node size |
orgDb |
perform clusterProfiler::setReadable based on this organism database |
bypassConverting |
bypass the symbol converting ID of rownames and those listed in EA result must be same |
edgeLink |
use geom_edge_link() instead of geom_edge_diagonal() |
cellLineName |
the cell line name to be included |
fontFamily |
font family name to be used for plotting |
strengthPlot |
append the barplot depicting edges with high strength |
nStrength |
specify how many edges are included in the strength plot |
strThresh |
the threshold for strength |
hub |
visualize the genes with top-n hub scores |
glowEdgeNum |
edges with top-n confidence of direction are highlighted |
nodePal |
vector of coloring of nodes (low, high) |
edgePal |
vector of coloring of edges (low, high) |
textCol |
color of texts in network plot |
titleCol |
color of title in network plot |
backCol |
color of background in network plot |
barTextCol |
text color in barplot |
barPal |
bar color |
barBackCol |
background color in barplot |
scoreType |
score type to use on inference |
barLegendKeyCol |
legend key color in barplot |
barAxisCol |
axis color in barplot |
bg.colour |
parameter to pass to geom_node_text |
bg.r |
parameter to pass to geom_node_text |
barPanelGridCol |
panel grid color in barplot |
titleSize |
the size of title |
seed |
A random seed to make the analysis reproducible, default is 1. |
ggplot2 object
data("exampleEaRes");data("exampleGeneExp") res <- bngeneplotCustom(results=exampleEaRes, exp=exampleGeneExp, pathNum=1, glowEdgeNum=NULL, hub=3, R=40, fontFamily="sans")
data("exampleEaRes");data("exampleGeneExp") res <- bngeneplotCustom(results=exampleEaRes, exp=exampleGeneExp, pathNum=1, glowEdgeNum=NULL, hub=3, R=40, fontFamily="sans")
Testing various R for bayesian network between genes
bngenetest( results, exp, expSample = NULL, algo = "hc", Rrange = seq(2, 40, 2), cl = NULL, algorithm.args = NULL, pathNum = NULL, convertSymbol = TRUE, expRow = "ENSEMBL", scoreType = "aic-g", orgDb = org.Hs.eg.db, bypassConverting = FALSE )
bngenetest( results, exp, expSample = NULL, algo = "hc", Rrange = seq(2, 40, 2), cl = NULL, algorithm.args = NULL, pathNum = NULL, convertSymbol = TRUE, expRow = "ENSEMBL", scoreType = "aic-g", orgDb = org.Hs.eg.db, bypassConverting = FALSE )
results |
the enrichment analysis result |
exp |
gene expression matrix |
expSample |
candidate rows to be included in the inference default to all |
algo |
structure learning method used in boot.strength() default to "hc" |
Rrange |
the sequence of R values to be tested |
cl |
cluster object from parallel::makeCluster() |
algorithm.args |
parameters to pass to bnlearn structure learnng function |
pathNum |
the pathway number (the number of row of the original result, ordered by p-value) |
convertSymbol |
whether the label of resulting network is converted to symbol, default to TRUE |
expRow |
the type of the identifier of rows of expression matrix |
scoreType |
return the specified scores |
orgDb |
perform clusterProfiler::setReadable based on this organism database |
bypassConverting |
bypass symbol converting |
list of graphs and scores
data("exampleEaRes");data("exampleGeneExp") res <- bngenetest(results = exampleEaRes, exp = exampleGeneExp, algo="hc", Rrange=seq(10, 30, 10), pathNum=1, scoreType="bge")
data("exampleEaRes");data("exampleGeneExp") res <- bngenetest(results = exampleEaRes, exp = exampleGeneExp, algo="hc", Rrange=seq(10, 30, 10), pathNum=1, scoreType="bge")
Plot pathway relationship
bnpathplot( results, exp, expSample = NULL, algo = "hc", algorithm.args = NULL, expRow = "ENSEMBL", cl = NULL, returnNet = FALSE, otherVar = NULL, otherVarName = NULL, qvalueCutOff = NULL, adjpCutOff = 0.05, nCategory = 15, R = 20, interactive = FALSE, color = "p.adjust", cexCategory = 1, cexLine = 0.5, chooseDir = FALSE, showDir = FALSE, delZeroDegree = TRUE, labelSize = 4, layout = "nicely", onlyDf = FALSE, disc = FALSE, tr = NULL, remainCont = NULL, shadowText = TRUE, bgColor = "white", textColor = "black", compareRef = FALSE, strThresh = NULL, strType = "normal", hub = NULL, scoreType = "bic-g", databasePal = "Set2", dep = NULL, sizeDep = FALSE, orgDb = org.Hs.eg.db, bypassConverting = FALSE, useSiGN = FALSE, edgeLink = TRUE, cellLineName = "5637_URINARY_TRACT", strengthPlot = FALSE, nStrength = 10, seed = 1 )
bnpathplot( results, exp, expSample = NULL, algo = "hc", algorithm.args = NULL, expRow = "ENSEMBL", cl = NULL, returnNet = FALSE, otherVar = NULL, otherVarName = NULL, qvalueCutOff = NULL, adjpCutOff = 0.05, nCategory = 15, R = 20, interactive = FALSE, color = "p.adjust", cexCategory = 1, cexLine = 0.5, chooseDir = FALSE, showDir = FALSE, delZeroDegree = TRUE, labelSize = 4, layout = "nicely", onlyDf = FALSE, disc = FALSE, tr = NULL, remainCont = NULL, shadowText = TRUE, bgColor = "white", textColor = "black", compareRef = FALSE, strThresh = NULL, strType = "normal", hub = NULL, scoreType = "bic-g", databasePal = "Set2", dep = NULL, sizeDep = FALSE, orgDb = org.Hs.eg.db, bypassConverting = FALSE, useSiGN = FALSE, edgeLink = TRUE, cellLineName = "5637_URINARY_TRACT", strengthPlot = FALSE, nStrength = 10, seed = 1 )
results |
the enrichment analysis result |
exp |
gene expression matrix |
expSample |
candidate rows to be included in the inference default to all |
algo |
structure learning method used in boot.strength() default to "hc" |
algorithm.args |
parameters to pass to bnlearn structure learnng function |
expRow |
the type of the identifier of rows of expression matrix |
cl |
cluster object from parallel::makeCluster() |
returnNet |
whether to return the network |
otherVar |
other variables to be included in the inference |
otherVarName |
the names of other variables |
qvalueCutOff |
the cutoff value for qvalue |
adjpCutOff |
the cutoff value for adjusted pvalues |
nCategory |
the number of pathways to be included |
R |
the number of bootstrap |
interactive |
whether to use bnviewer (default to FALSE) |
color |
color of node, default to adjusted p-value |
cexCategory |
scaling factor of size of nodes |
cexLine |
scaling factor of width of edges |
chooseDir |
if undirected edges are present, choose direction of edges |
showDir |
show the confidence of direction of edges |
delZeroDegree |
delete zero degree nodes |
labelSize |
the size of label of the nodes |
layout |
ggraph layout, default to "nicely" |
onlyDf |
return only data.frame used for inference |
disc |
discretize the expressoin data |
tr |
Specify data.frame if one needs to discretize as the same parameters as the other dataset |
remainCont |
Specify characters when perform discretization, if some columns are to be remain continuous |
shadowText |
whether to use shadow text for the better readability (default: TRUE) |
bgColor |
color for text background when shadowText is TRUE |
textColor |
color for text when shadowText is TRUE |
compareRef |
whether compare to the reference network between pathway |
strThresh |
threshold for strength, automatically determined if NULL |
strType |
"normal" or "ms" for multiscale implementation |
hub |
change the shape of node according to hub scores (default NULL) |
scoreType |
score type to use on choosing edge direction |
databasePal |
palette to be used in scale_color_brewer when the multiple results are to be shown |
dep |
the tibble storing dependency score from library depmap |
sizeDep |
whether to reflect DepMap score to the node size |
orgDb |
perform clusterProfiler::setReadable based on this organism database |
bypassConverting |
bypass the symbol converting If you use custom annotation databases that does not have SYMBOL listed in keys. ID of rownames and those listed in EA result must be same. |
useSiGN |
default to FALSE. For using SiGN-BN in the function in Windows 10/11, 1. Download the SiGN-BN HC+BS binary in WSL (https://sign.hgc.jp/signbn/download.html) 2. Set PATH to executable (sign.1.8.3) |
edgeLink |
whether to set edge to geom_edge_link() FALSE to use geom_edge_diagonal() |
cellLineName |
the cell line name to be included |
strengthPlot |
append the barplot depicting edges with high strength |
nStrength |
specify how many edges are included in the strength plot |
seed |
A random seed to make the analysis reproducible, default is 1. |
ggplot2 object
data("exampleEaRes");data("exampleGeneExp") res <- bnpathplot(results = exampleEaRes, exp = exampleGeneExp, R = 10, expRow = "ENSEMBL")
data("exampleEaRes");data("exampleGeneExp") res <- bnpathplot(results = exampleEaRes, exp = exampleGeneExp, R = 10, expRow = "ENSEMBL")
Plot pathway relationship using customized theme
bnpathplotCustom( results, exp, expSample = NULL, algo = "hc", R = 20, expRow = "ENSEMBL", color = "p.adjust", cexCategory = 1, cl = NULL, showDir = FALSE, chooseDir = FALSE, labelSize = 4, layout = "nicely", strType = "normal", compareRef = FALSE, disc = FALSE, tr = NULL, remainCont = NULL, qvalueCutOff = NULL, adjpCutOff = 0.05, nCategory = 15, cexLine = 1, returnNet = FALSE, dep = NULL, sizeDep = FALSE, cellLineName = "5637_URINARY_TRACT", fontFamily = "sans", otherVar = NULL, otherVarName = NULL, onlyDf = FALSE, algorithm.args = NULL, strengthPlot = FALSE, nStrength = 10, edgeLink = FALSE, strThresh = NULL, hub = NULL, glowEdgeNum = NULL, nodePal = c("blue", "red"), edgePal = c("blue", "red"), textCol = "black", backCol = "white", barTextCol = "black", barPal = c("red", "blue"), barBackCol = "white", scoreType = "bic-g", barLegendKeyCol = "white", orgDb = org.Hs.eg.db, barAxisCol = "black", barPanelGridCol = "black", bg.colour = NULL, bg.r = 0.1, seed = 1, bypassConverting = FALSE )
bnpathplotCustom( results, exp, expSample = NULL, algo = "hc", R = 20, expRow = "ENSEMBL", color = "p.adjust", cexCategory = 1, cl = NULL, showDir = FALSE, chooseDir = FALSE, labelSize = 4, layout = "nicely", strType = "normal", compareRef = FALSE, disc = FALSE, tr = NULL, remainCont = NULL, qvalueCutOff = NULL, adjpCutOff = 0.05, nCategory = 15, cexLine = 1, returnNet = FALSE, dep = NULL, sizeDep = FALSE, cellLineName = "5637_URINARY_TRACT", fontFamily = "sans", otherVar = NULL, otherVarName = NULL, onlyDf = FALSE, algorithm.args = NULL, strengthPlot = FALSE, nStrength = 10, edgeLink = FALSE, strThresh = NULL, hub = NULL, glowEdgeNum = NULL, nodePal = c("blue", "red"), edgePal = c("blue", "red"), textCol = "black", backCol = "white", barTextCol = "black", barPal = c("red", "blue"), barBackCol = "white", scoreType = "bic-g", barLegendKeyCol = "white", orgDb = org.Hs.eg.db, barAxisCol = "black", barPanelGridCol = "black", bg.colour = NULL, bg.r = 0.1, seed = 1, bypassConverting = FALSE )
results |
the enrichment analysis result |
exp |
gene expression matrix |
expSample |
candidate rows to be included in the inference default to all |
algo |
structure learning method used in boot.strength() default to "hc" |
R |
the number of bootstrap |
expRow |
the type of the identifier of rows of expression matrix |
color |
color of node, default to adjusted p-value |
cexCategory |
scaling factor of size of nodes |
cl |
cluster object from parallel::makeCluster() |
showDir |
show the confidence of direction of edges |
chooseDir |
if undirected edges are present, choose direction of edges |
labelSize |
the size of label of the nodes |
layout |
ggraph layout, default to "nicely" |
strType |
"normal" or "ms" for multiscale implementation |
compareRef |
whether compare to the reference network between pathway |
disc |
discretize the expressoin data |
tr |
Specify data.frame if one needs to discretize as the same parameters as the other dataset |
remainCont |
Specify characters when perform discretization, if some columns are to be remain continuous |
qvalueCutOff |
the cutoff value for qvalue |
adjpCutOff |
the cutoff value for adjusted pvalues |
nCategory |
the number of pathways to be included |
cexLine |
scaling factor of width of edges |
returnNet |
whether to return the network |
dep |
the tibble storing dependency score from library depmap |
sizeDep |
whether to reflect DepMap score to the node size |
cellLineName |
the cell line name to be included |
fontFamily |
font family name to be used for plotting |
otherVar |
other variables to be included in the inference |
otherVarName |
the names of other variables |
onlyDf |
return only data.frame used for inference |
algorithm.args |
parameters to pass to bnlearn structure learnng function |
strengthPlot |
append the barplot depicting edges with high strength |
nStrength |
specify how many edges are included in the strength plot |
edgeLink |
use geom_edge_link() instead of geom_edge_diagonal() |
strThresh |
threshold for strength, automatically determined if NULL |
hub |
change the shape of node according to hub scores (default NULL) |
glowEdgeNum |
edges with top-n confidence of direction are highlighted |
nodePal |
vector of coloring of nodes (low, high) |
edgePal |
vector of coloring of edges (low, high) |
textCol |
color of texts in network plot |
backCol |
color of background in network plot |
barTextCol |
text color in barplot |
barPal |
bar color |
barBackCol |
background color in barplot |
scoreType |
score type to use on inference |
barLegendKeyCol |
legend key color in barplot |
orgDb |
perform clusterProfiler::setReadable based on this organism database |
barAxisCol |
axis color in barplot |
barPanelGridCol |
panel grid color in barplot |
bg.colour |
parameter to pass to geom_node_text |
bg.r |
parameter to pass to geom_node_text |
seed |
A random seed to make the analysis reproducible, default is 1. |
bypassConverting |
bypass the symbol converting ID of rownames and those listed in EA result must be same |
ggplot2 object
data("exampleEaRes");data("exampleGeneExp") res <- bnpathplotCustom(results=exampleEaRes, exp=exampleGeneExp, fontFamily="sans", glowEdgeNum=3, hub=3)
data("exampleEaRes");data("exampleGeneExp") res <- bnpathplotCustom(results=exampleEaRes, exp=exampleGeneExp, fontFamily="sans", glowEdgeNum=3, hub=3)
Testing various R for bayesian network between pathways
bnpathtest( results, exp, expSample = NULL, algo = "hc", algorithm.args = NULL, expRow = "ENSEMBL", cl = NULL, orgDb = org.Hs.eg.db, bypassConverting = FALSE, qvalueCutOff = 0.05, adjpCutOff = 0.05, nCategory = 15, Rrange = seq(2, 40, 2), scoreType = "aic-g" )
bnpathtest( results, exp, expSample = NULL, algo = "hc", algorithm.args = NULL, expRow = "ENSEMBL", cl = NULL, orgDb = org.Hs.eg.db, bypassConverting = FALSE, qvalueCutOff = 0.05, adjpCutOff = 0.05, nCategory = 15, Rrange = seq(2, 40, 2), scoreType = "aic-g" )
results |
the enrichment analysis result |
exp |
gene expression matrix |
expSample |
candidate rows to be included in the inference default to all |
algo |
structure learning method used in boot.strength() default to "hc" |
algorithm.args |
parameters to pass to bnlearn structure learnng function |
expRow |
the type of the identifier of rows of expression matrix |
cl |
cluster object from parallel::makeCluster() |
orgDb |
perform clusterProfiler::setReadable based on this organism database |
bypassConverting |
bypass symbol converting |
qvalueCutOff |
the cutoff value for qvalue |
adjpCutOff |
the cutoff value for adjusted pvalues |
nCategory |
the number of pathways to be included |
Rrange |
the sequence of R values to be tested |
scoreType |
return the specified scores |
list of graphs and scores
data("exampleEaRes");data("exampleGeneExp") res <- bnpathtest(results = exampleEaRes, exp = exampleGeneExp, algo="hc", Rrange=seq(10, 30, 10), expRow = "ENSEMBL", scoreType="bge")
data("exampleEaRes");data("exampleGeneExp") res <- bnpathtest(results = exampleEaRes, exp = exampleGeneExp, algo="hc", Rrange=seq(10, 30, 10), expRow = "ENSEMBL", scoreType="bge")
Take the list of networks and returns the F-measures
compareBNs(listOfNets)
compareBNs(listOfNets)
listOfNets |
list of networks |
F-measures of each combination of network
data("exampleEaRes");data("exampleGeneExp") net1 <- bngeneplot(results = exampleEaRes, exp = exampleGeneExp, pathNum = 1, R = 10, returnNet=TRUE) net2 <- bngeneplot(results = exampleEaRes, exp = exampleGeneExp, pathNum = 1, R = 10, returnNet=TRUE) res <- compareBNs(list(net1$av, net2$av))
data("exampleEaRes");data("exampleGeneExp") net1 <- bngeneplot(results = exampleEaRes, exp = exampleGeneExp, pathNum = 1, R = 10, returnNet=TRUE) net2 <- bngeneplot(results = exampleEaRes, exp = exampleGeneExp, pathNum = 1, R = 10, returnNet=TRUE) res <- compareBNs(list(net1$av, net2$av))
An example enrichment analysis result to be used for testing purpose. The result was produced by running ReactomePA::enrichPathway() and subsequent clusterProfiler::setReadable() on 'exampleGeneExp'.
data(exampleEaRes)
data(exampleEaRes)
An object of class enrichResult
with 47 rows and 9 columns.
example enrichment analysis result
An example gene expression data to be used for testing purpose made by runif() for ERCC genes and 100 samples. No biological meanings can be obtained from the data.
data(exampleGeneExp)
data(exampleGeneExp)
An object of class data.frame
with 7 rows and 100 columns.
example gene expression
multiscale bootstrap-based inference of Bayesian network
inferMS(data, algo, algorithm.args, R, cl = NULL, r = seq(0.5, 1.5, 0.1))
inferMS(data, algo, algorithm.args, R, cl = NULL, r = seq(0.5, 1.5, 0.1))
data |
data.frame to perform inference |
algo |
structure learning method used in boot.strength() |
algorithm.args |
parameters to pass to bnlearn structure learnng function |
R |
the number of bootstrap |
cl |
cluster object from parallel::makeCluster() |
r |
vector for size of each bootstrap replicate |
object of class bn.strength
Load the output of SiGN-BN (HC+BS)
loadSign(fileName)
loadSign(fileName)
fileName |
the result of SiGN-BN |
list of edges, nodes, strength, and bn (bnlearn)
obtain the analysis results including the queried gene symbol
obtainPath(res, geneSymbol)
obtainPath(res, geneSymbol)
res |
enrichment analysis result |
geneSymbol |
the candidate gene |
subset of enrichment results
data("exampleEaRes") obtainPath(res = exampleEaRes, geneSymbol="ERCC7")
data("exampleEaRes") obtainPath(res = exampleEaRes, geneSymbol="ERCC7")
produce a plot of bnlearn::cpdist by performing bnlearn::cpdist on specified node, evidence and level.
queryCpDistLs(fitted, candidate, evidences, discPalette = "Set2", ...)
queryCpDistLs(fitted, candidate, evidences, discPalette = "Set2", ...)
fitted |
bn.fit object |
candidate |
name of node |
evidences |
the evidences |
discPalette |
palette to be used for plotting if the event is discrete |
... |
other parameters passed to bnlearn cpdist |
list of dataframe containing raw values
library(bnlearn) data("exampleEaRes") data("exampleGeneExp") net <- bngeneplot(exampleEaRes, exampleGeneExp, pathNum=1, returnNet=TRUE) fitted <- bn.fit(net$av, net$df) res <- queryCpDistLs(fitted, candidate="ERCC4", evidences=c("ERCC2<0.1","ERCC2>0.5","ERCC2>0.8"), n=500)
library(bnlearn) data("exampleEaRes") data("exampleGeneExp") net <- bngeneplot(exampleEaRes, exampleGeneExp, pathNum=1, returnNet=TRUE) fitted <- bn.fit(net$av, net$df) res <- queryCpDistLs(fitted, candidate="ERCC4", evidences=c("ERCC2<0.1","ERCC2>0.5","ERCC2>0.8"), n=500)
produce a plot of bnlearn::cpdist by performing bnlearn::cpdist on specified node, evidence and level.
queryCpDistLw( fitted, candidate, evidence, levels, point = FALSE, pointSize = 5, alpha = TRUE, ... )
queryCpDistLw( fitted, candidate, evidence, levels, point = FALSE, pointSize = 5, alpha = TRUE, ... )
fitted |
bn.fit object |
candidate |
name of node |
evidence |
evidence variable name |
levels |
level to be listed |
point |
geom_point the weighted mean |
pointSize |
point size for geom_point |
alpha |
whether to reflect the weights by alpha (TRUE) or color (FALSE) |
... |
other parameters passed to bnlearn cpdist |
list of dataframe containing raw values
library(bnlearn) data("exampleEaRes") data("exampleGeneExp") net <- bngeneplot(exampleEaRes, exampleGeneExp, pathNum=1, returnNet=TRUE) fitted <- bn.fit(net$av, net$df) res <- queryCpDistLw(fitted, candidate="ERCC4", evidence="ERCC2", levels=c(0.1, 0.5, 0.8), n=500)
library(bnlearn) data("exampleEaRes") data("exampleGeneExp") net <- bngeneplot(exampleEaRes, exampleGeneExp, pathNum=1, returnNet=TRUE) fitted <- bn.fit(net$av, net$df) res <- queryCpDistLw(fitted, candidate="ERCC4", evidence="ERCC2", levels=c(0.1, 0.5, 0.8), n=500)