Title: | Combines WGCNA and xCell readouts with bayesian network learrning to generate a Gene-Module Immune-Cell network (GMIC) |
---|---|
Description: | This package uses bayesian network learning to detect relationships between Gene Modules detected by WGCNA and immune cell signatures defined by xCell. It is a hypothesis generating tool. |
Authors: | Richard Virgen-Slane |
Maintainer: | Richard Virgen-Slane <[email protected]> |
License: | GPL-2 + file LICENSE |
Version: | 1.21.0 |
Built: | 2024-10-30 09:16:53 UTC |
Source: | https://github.com/bioc/GmicR |
Carries out WGCNA with default settings or custom settings
Auto_WGCNA( datExpr, colname_correct = TRUE, minModuleSize = 10, deepSplit = 4, networkType = "signed hybrid", TOMType = "unsigned", corFnc = "bicor", mergeCutHeight = 0.25, sft_RsquaredCut = 0.85, removeFirst = FALSE, reassignThreshold = 1e-06, maxBlockSize = 25000, nThreads = NULL )
Auto_WGCNA( datExpr, colname_correct = TRUE, minModuleSize = 10, deepSplit = 4, networkType = "signed hybrid", TOMType = "unsigned", corFnc = "bicor", mergeCutHeight = 0.25, sft_RsquaredCut = 0.85, removeFirst = FALSE, reassignThreshold = 1e-06, maxBlockSize = 25000, nThreads = NULL )
datExpr |
Expression data. A matrix (preferred) or
data frame in which columns are genes and rows ar samples. NAs are
allowed, but not too many. See |
colname_correct |
a logical value. If TRUE (default), "." in gene names will be replaced with "-". This corrects a name change that is induced by R when creating a data.frame. If FALSE, no changes will be made. |
minModuleSize |
minimum module size for module detection. See
|
deepSplit |
integer value between 0 and 4. Provides a simplified control over how sensitive
module detection should be to module splitting, with 0 least and 4 most sensitive. See
|
networkType |
network type. Allowed values are (unique abbreviations of) |
TOMType |
one of |
corFnc |
the correlation function to be used in adjacency calculation. |
mergeCutHeight |
dendrogram cut height for module merging. |
sft_RsquaredCut |
desired minimum scale free topology fitting index R^2. Default is 0.80. |
removeFirst |
should the first bin be removed from the connectivity histogram? |
reassignThreshold |
p-value ratio threshold for reassigning genes between modules. See Details. |
maxBlockSize |
integer giving maximum block size for module detection. Ignored if |
nThreads |
non-negative integer specifying the number of parallel threads to be used by certain parts of correlation calculations. This option only has an effect on systems on which a POSIX thread library is available (which currently includes Linux and Mac OSX, but excludes Windows). If zero, the number of online processors will be used if it can be determined dynamically, otherwise correlation calculations will use 2 threads. |
Returns a lists containing network input parameters used for WGCNA, WGCNA module information, and quality control plots.
This is a wrapper for WGCNA.
sample_dat_dir<-system.file("extdata", "sample_dat.Rdata", package = "GmicR", mustWork = TRUE) load(sample_dat_dir) GMIC_Builder<-Auto_WGCNA(sample_dat, mergeCutHeight = 0.35, minModuleSize = 10)
sample_dat_dir<-system.file("extdata", "sample_dat.Rdata", package = "GmicR", mustWork = TRUE) load(sample_dat_dir) GMIC_Builder<-Auto_WGCNA(sample_dat, mergeCutHeight = 0.35, minModuleSize = 10)
Generates a subgraph from query nodes
Batch_Net(bn_output, Node_ids, relationship_type = "nbr")
Batch_Net(bn_output, Node_ids, relationship_type = "nbr")
bn_output |
R object output from bn_tabu_gen function |
Node_ids |
vector containing the nodes for subgraph generation. |
relationship_type |
the relationship to be returned for the specified query nodes. The options are "mb", "nbr", "parents", "children". Default setting is "nbr". |
a subgraph containg the selected nodes and relationships.
Uses tabu search algorithm to learn the structure of discretized data.
bn_tabu_gen( Auto_WGCNA_OUTPUT, whitelist = NULL, blacklist = NULL, score = "bde", tabu = 50, iss = 10, cluster = NULL, debug = TRUE, bootstraps_replicates = 500 )
bn_tabu_gen( Auto_WGCNA_OUTPUT, whitelist = NULL, blacklist = NULL, score = "bde", tabu = 50, iss = 10, cluster = NULL, debug = TRUE, bootstraps_replicates = 500 )
Auto_WGCNA_OUTPUT |
an R object generated by Auto_WGCNA and discretized using the Data_Prep function. |
whitelist |
a data frame with two columns (optionally labeled "from" and "to"), containing a set of arcs to be included in the graph. |
blacklist |
a data frame with two columns (optionally labeled "from" and "to"), containing a set of arcs not to be included in the graph. |
score |
character string indicating the score used for structure learning. If "bde" (default), prior is set to "uniform". If bds is used, the prior is set to "marginal". |
tabu |
a positive integer number, the length of the tabu list used in the
|
iss |
the imaginary sample size, used by the Bayesian Dirichlet scores (bde and bds) It is also known as “equivalent sample size”. The default value is equal to 10. |
cluster |
an optional cluster object from package parallel. |
debug |
a boolean value. If |
bootstraps_replicates |
an integer for the number of bootstraps_replicates used for structure learning. Default value is 500 |
The learned bayesian network
GMIC_Builder_disc_dir<-system.file("extdata", "GMIC_Builder_disc.Rdata", package = "GmicR", mustWork = TRUE) load(GMIC_Builder_disc_dir) no_cores<-1 cl<-parallel::makeCluster(no_cores) GMIC_net<-bn_tabu_gen(GMIC_Builder_disc, cluster = cl, bootstraps_replicates = 50, score = "bds") parallel::stopCluster(cl)
GMIC_Builder_disc_dir<-system.file("extdata", "GMIC_Builder_disc.Rdata", package = "GmicR", mustWork = TRUE) load(GMIC_Builder_disc_dir) no_cores<-1 cl<-parallel::makeCluster(no_cores) GMIC_net<-bn_tabu_gen(GMIC_Builder_disc, cluster = cl, bootstraps_replicates = 50, score = "bds") parallel::stopCluster(cl)
Discretizes biological assay data in preparation for bayensian network learning
Data_Prep( Auto_WGCNA_OUTPUT = NULL, Remove_ME0 = FALSE, Numeric_Pheno_scores = NULL, xCell_Signatures = NULL, ibreaks = 60 )
Data_Prep( Auto_WGCNA_OUTPUT = NULL, Remove_ME0 = FALSE, Numeric_Pheno_scores = NULL, xCell_Signatures = NULL, ibreaks = 60 )
Auto_WGCNA_OUTPUT |
R object generated from Auto_WGCNA function. |
Remove_ME0 |
a logical value. If FALSE (default), ME0 is not removed. If TRUE the eigengene for module 0 is removed prior to analysis. |
Numeric_Pheno_scores |
a data.frame with rows indicating sample ID and columns representing additional phenotype data to be included in BN learning. If NULL (default) no data will be included. If provided, the data.frame will be merged with MEs and discretized into three levels. |
xCell_Signatures |
the name of the text file generated by xCell that contains the cell signature scores. If NULL (default) the only module eigenegnes will be processed. If not NULL and if Auto_WGCNA_OUTPUT is NULL, cell signature scores will be discretized. |
ibreaks |
an integer that indicates the number of ibreaks used for discretization. The default value is 60. |
a list containing a data.frame with module eigenegnes merged with Xcell signature scores and discretized into three levels: L, M, H. If Auto_WGCNA_OUTPUT is NULL, both scaled and discretized cell signatures will be return.
Please verify that the sample name formatting is consistent between both datasets. Rownames in the module eigengenes data.frame and the column names of xCell signatures scores text file are matched for merging. Only samples that are present in both will be processed!
file_dir<-system.file("extdata", "IRIS_xCell_sig.txt", package = "GmicR", mustWork = TRUE) Disc_Xcell_sig<-Data_Prep(xCell_Signatures=file_dir, ibreaks = 10) Disc_Xcell_sig$disc_data
file_dir<-system.file("extdata", "IRIS_xCell_sig.txt", package = "GmicR", mustWork = TRUE) Disc_Xcell_sig<-Data_Prep(xCell_Signatures=file_dir, ibreaks = 10) Disc_Xcell_sig$disc_data
Visualized network
Gmic_viz(Auto_WGCNA_Output, Filter_unconnected_ME = TRUE)
Gmic_viz(Auto_WGCNA_Output, Filter_unconnected_ME = TRUE)
Auto_WGCNA_Output |
R object with GMIC bayesian network |
Filter_unconnected_ME |
a logical value. If TRUE, the default, unconnected modules will be removed from the final network. If FALSE, all modules will be shown. |
a shiny object for network visualization.
GMIC_Final_dir<-system.file("extdata", "GMIC_Final.Rdata", package = "GmicR", mustWork = TRUE) load(GMIC_Final_dir) if(interactive()){ Gmic_viz(GMIC_Final)}
GMIC_Final_dir<-system.file("extdata", "GMIC_Final.Rdata", package = "GmicR", mustWork = TRUE) load(GMIC_Final_dir) if(interactive()){ Gmic_viz(GMIC_Final)}
GO enrichment for module names
GO_Module_NameR(Auto_WGCNA_OUTPUT, cutoff_size = NULL)
GO_Module_NameR(Auto_WGCNA_OUTPUT, cutoff_size = NULL)
Auto_WGCNA_OUTPUT |
object returned by GSEAGO_Builder function |
cutoff_size |
integer indication the GO size cut off. If NULL (default) the term with the smallest Pvalue will be returned. |
a data.frame listing the GO name for each module. Additionally, a table similar to the output of Query_Prep function is returned, appended with module GO names.
GMIC_Builder_dir<-system.file("extdata", "GMIC_Builder.Rdata", package = "GmicR", mustWork = TRUE) load(GMIC_Builder_dir) GMIC_Builder<-GO_Module_NameR(GMIC_Builder, cutoff_size = 100)
GMIC_Builder_dir<-system.file("extdata", "GMIC_Builder.Rdata", package = "GmicR", mustWork = TRUE) load(GMIC_Builder_dir) GMIC_Builder<-GO_Module_NameR(GMIC_Builder, cutoff_size = 100)
GO enrichment for module names
GSEAGO_Builder( Auto_WGCNA_OUTPUT, species = "Homo sapiens", no_cores = 4, ontology = "BP", GO_conditional = FALSE, colname_correct = TRUE )
GSEAGO_Builder( Auto_WGCNA_OUTPUT, species = "Homo sapiens", no_cores = 4, ontology = "BP", GO_conditional = FALSE, colname_correct = TRUE )
Auto_WGCNA_OUTPUT |
output from Auto_WGCNA function. |
species |
either 'Homo sapiens' (default) or 'Mus musculus'. |
no_cores |
Number of cores to use. Default = 4. |
ontology |
string either 'BP'(Biological Process; default), 'CC'(Cellular Component), or 'MF' (Molecular Function). |
GO_conditional |
A logical indicating whether the calculation should condition on the GO structure. will not be carried out. If TRUE, |
colname_correct |
a logical value. If TRUE (default), "." in gene names will be replaced with "-". This corrects a name change that is induced by R when creating a data.frame. If FALSE, no changes will be made. |
Lists with gene ontology enrichment analysis, performed using GOstats, for each module.
gene names must be official gene symbol
GMIC_Builder_dir<-system.file("extdata", "GMIC_Builder.Rdata", package = "GmicR", mustWork = TRUE) load(GMIC_Builder_dir) GMIC_Builder$GSEAGO_Builder_Output<-NULL Test_GMIC_Builder<-GSEAGO_Builder(GMIC_Builder, no_cores = 1) summary(Test_GMIC_Builder$GSEAGO_Builder_Output)
GMIC_Builder_dir<-system.file("extdata", "GMIC_Builder.Rdata", package = "GmicR", mustWork = TRUE) load(GMIC_Builder_dir) GMIC_Builder$GSEAGO_Builder_Output<-NULL Test_GMIC_Builder<-GSEAGO_Builder(GMIC_Builder, no_cores = 1) summary(Test_GMIC_Builder$GSEAGO_Builder_Output)
Identifies arcs between nodes with inverse relationships
InverseARCs(Output, threshold = -0.3)
InverseARCs(Output, threshold = -0.3)
Output |
a data frame containing the output of BN_Conditions function. |
threshold |
number indicating the maximum slope for defining negative relationships. Default level is -0.3. |
arcs with inverse relationships
GMIC_net_dir<-system.file("extdata", "GMIC_net.Rdata", package = "GmicR", mustWork = TRUE) load(GMIC_net_dir) GMIC_Final<-InverseARCs(GMIC_net, threshold = -0.3)
GMIC_net_dir<-system.file("extdata", "GMIC_net.Rdata", package = "GmicR", mustWork = TRUE) load(GMIC_net_dir) GMIC_Final<-InverseARCs(GMIC_net, threshold = -0.3)
Query Prep
Query_Prep( Auto_WGCNA_OUTPUT, numGenes = 500, Find_hubs = FALSE, nThreads = NULL, calculate_intramodularConnectivity = TRUE )
Query_Prep( Auto_WGCNA_OUTPUT, numGenes = 500, Find_hubs = FALSE, nThreads = NULL, calculate_intramodularConnectivity = TRUE )
Auto_WGCNA_OUTPUT |
R object generated by Auto_WGCNA function. |
numGenes |
integer indicating the number of random genes to test for hub gene detection. Default is 500. |
Find_hubs |
logical value. If TRUE, module hub genes will be returned. If FALSE (default), intramodularConnectivity will be returned without hub gene identification. |
nThreads |
non-negative integer specifying the number of parallel threads to be used by certain parts of correlation calculations. This option only has an effect on systems on which a POSIX thread library is available (which currently includes Linux and Mac OSX, but excludes Windows). If zero, the number of online processors will be used if it can be determined dynamically, otherwise correlation calculations will use 2 threads. |
calculate_intramodularConnectivity |
a logical value. If TRUE (default), the intramodularConnectivity will be caluculated using the intramodularConnectivity function from WGCNA. If FALSE, a table of modules and genes will be returned without intramodularConnectivity values. |
a data.frame detailing the gene symbols for each module. Gene intramodularConnectivity may also be returned. If detected, hub genes are annotated.
GMIC_Builder_dir<-system.file("extdata", "GMIC_Builder.Rdata", package = "GmicR", mustWork = TRUE) load(GMIC_Builder_dir) GMIC_Builder<-Query_Prep(GMIC_Builder, Find_hubs = TRUE) head(GMIC_Builder$Query)
GMIC_Builder_dir<-system.file("extdata", "GMIC_Builder.Rdata", package = "GmicR", mustWork = TRUE) load(GMIC_Builder_dir) GMIC_Builder<-Query_Prep(GMIC_Builder, Find_hubs = TRUE) head(GMIC_Builder$Query)
Scales and centers data by sample/row in preparation for discretization
xCell_loader(File = NULL)
xCell_loader(File = NULL)
File |
the name of the text file generated by xCell that contains the cell signature scores. |
xCell signatures scaled and centered by sample. For GMIC, ImmuneScore, StromaScore, and MicroenvironmentScore are removed.
file_dir<-system.file("extdata", "IRIS_xCell_sig.txt", package = "GmicR", mustWork = TRUE) Xcell_sig<-xCell_loader(file_dir) plot(Xcell_sig$Bcells)
file_dir<-system.file("extdata", "IRIS_xCell_sig.txt", package = "GmicR", mustWork = TRUE) Xcell_sig<-xCell_loader(file_dir) plot(Xcell_sig$Bcells)