| trans_network {microeco} | R Documentation | 
Create trans_network object for network analysis.
Description
This class is a wrapper for a series of network analysis methods, including the network construction, network attributes analysis, eigengene analysis, network subsetting, node and edge properties, network visualization and other operations.
Methods
Public methods
Method new()
Create the trans_network object, store the important intermediate data 
and calculate correlations if cor_method parameter is not NULL.
Usage
trans_network$new(
  dataset = NULL,
  cor_method = NULL,
  use_WGCNA_pearson_spearman = FALSE,
  use_NetCoMi_pearson_spearman = FALSE,
  use_sparcc_method = c("NetCoMi", "SpiecEasi")[1],
  taxa_level = "OTU",
  filter_thres = 0,
  nThreads = 1,
  SparCC_simu_num = 100,
  env_cols = NULL,
  add_data = NULL,
  ...
)Arguments
- dataset
- default NULL; the object of - microtableclass. Default NULL means customized analysis.
- cor_method
- default NULL; NULL or one of "bray", "pearson", "spearman", "sparcc", "bicor", "cclasso" and "ccrepe"; All the methods refered to - NetCoMipackage are performed based on- netConstructfunction of- NetCoMipackage and require- NetCoMito be installed from Github (https://github.com/stefpeschel/NetCoMi); For the algorithm details, please see Peschel et al. 2020 Brief. Bioinform <doi: 10.1093/bib/bbaa290>;- NULL
- NULL denotes non-correlation network, i.e. do not use correlation-based network. If so, the return res_cor_p list will be NULL. 
- 'bray'
- 1-B, where B is Bray-Curtis dissimilarity; based on - vegan::vegdistfunction
- 'pearson'
- Pearson correlation; If - use_WGCNA_pearson_spearmanand- use_NetCoMi_pearson_spearmanare both FALSE, use the function- cor.testin R;- use_WGCNA_pearson_spearman = TRUEinvoke- corAndPvaluefunction of- WGCNApackage;- use_NetCoMi_pearson_spearman = TRUEinvoke- netConstructfunction of- NetCoMipackage
- 'spearman'
- Spearman correlation; other details are same with the 'pearson' option 
- 'sparcc'
- SparCC algorithm (Friedman & Alm, PLoS Comp Biol, 2012, <doi:10.1371/journal.pcbi.1002687>); use NetCoMi package when - use_sparcc_method = "NetCoMi"; use- SpiecEasipackage when- use_sparcc_method = "SpiecEasi"and require- SpiecEasito be installed from Github (https://github.com/zdk123/SpiecEasi)
- 'bicor'
- Calculate biweight midcorrelation efficiently for matrices based on - WGCNA::bicorfunction; This option can invoke- netConstructfunction of- NetCoMipackage; Make sure- WGCNAand- NetCoMipackages are both installed
- 'cclasso'
- Correlation inference of Composition data through Lasso method based on - netConstructfunction of- NetCoMipackage; for details, see- NetCoMi::cclassofunction
- 'ccrepe'
- Calculates compositionality-corrected p-values and q-values for compositional data using an arbitrary distance metric based on - NetCoMi::netConstructfunction; also see- NetCoMi::ccrepefunction
 
- use_WGCNA_pearson_spearman
- default FALSE; whether use WGCNA package to calculate correlation when - cor_method= "pearson" or "spearman".
- use_NetCoMi_pearson_spearman
- default FALSE; whether use NetCoMi package to calculate correlation when - cor_method= "pearson" or "spearman". The important difference between NetCoMi and others is the features of zero handling and data normalization; See <doi: 10.1093/bib/bbaa290>.
- use_sparcc_method
- default - c("NetCoMi", "SpiecEasi")[1]; use- NetCoMipackage or- SpiecEasipackage to perform SparCC when- cor_method = "sparcc".
- taxa_level
- default "OTU"; taxonomic rank; 'OTU' denotes using feature abundance table; other available options should be one of the colnames of - tax_tableof input dataset.
- filter_thres
- default 0; the relative abundance threshold. 
- nThreads
- default 1; the CPU thread number; available when - use_WGCNA_pearson_spearman = TRUEor- use_sparcc_method = "SpiecEasi".
- SparCC_simu_num
- default 100; SparCC simulation number for bootstrap when - use_sparcc_method = "SpiecEasi".
- env_cols
- default NULL; numeric or character vector to select the column names of environmental data in dataset$sample_table; the environmental data can be used in the correlation network (as the nodes) or - FlashWeavenetwork.
- add_data
- default NULL; provide environmental variable table additionally instead of - env_colsparameter; rownames must be sample names.
- ...
- parameters pass to - NetCoMi::netConstructfor other operations, such as zero handling and/or data normalization when cor_method and other parameters refer to- NetCoMipackage.
Returns
res_cor_p list with the correlation (association) matrix and p value matrix. Note that when cor_method and other parameters
refer to NetCoMi package, the p value table are all zero as the significant associations have been selected.
Examples
\donttest{
data(dataset)
# for correlation network
t1 <- trans_network$new(dataset = dataset, cor_method = "pearson", 
		taxa_level = "OTU", filter_thres = 0.0002)
# for non-correlation network
t1 <- trans_network$new(dataset = dataset, cor_method = NULL)
}
Method cal_network()
Construct network based on the igraph package or SpiecEasi package or julia FlashWeave package or beemStatic package.
Usage
trans_network$cal_network(
  network_method = c("COR", "SpiecEasi", "gcoda", "FlashWeave", "beemStatic")[1],
  COR_p_thres = 0.01,
  COR_p_adjust = "fdr",
  COR_weight = TRUE,
  COR_cut = 0.6,
  COR_optimization = FALSE,
  COR_optimization_low_high = c(0.01, 0.8),
  COR_optimization_seq = 0.01,
  SpiecEasi_method = "mb",
  FlashWeave_tempdir = NULL,
  FlashWeave_meta_data = FALSE,
  FlashWeave_other_para = "alpha=0.01,sensitive=true,heterogeneous=true",
  beemStatic_t_strength = 0.001,
  beemStatic_t_stab = 0.8,
  add_taxa_name = "Phylum",
  delete_unlinked_nodes = TRUE,
  usename_rawtaxa_when_taxalevel_notOTU = FALSE,
  ...
)Arguments
- network_method
- default "COR"; "COR", "SpiecEasi", "gcoda", "FlashWeave" or "beemStatic"; - network_method = NULLmeans skipping the network construction for the customized use. The option details:- 'COR'
- correlation-based network; use the correlation and p value matrices in - res_cor_plist stored in the object; See Deng et al. (2012) <doi:10.1186/1471-2105-13-113> for other details
- 'SpiecEasi'
- SpiecEasinetwork; relies on algorithms of sparse neighborhood and inverse covariance selection; belong to the category of conditional dependence and graphical models; see https://github.com/zdk123/SpiecEasi for installing the R package; see Kurtz et al. (2015) <doi:10.1371/journal.pcbi.1004226> for the algorithm details
- 'gcoda'
- hypothesize the logistic normal distribution of microbiome data; use penalized maximum likelihood method to estimate the sparse structure of inverse covariance for latent normal variables to address the high dimensionality of the microbiome data; belong to the category of conditional dependence and graphical models; depend on the R - NetCoMipackage https://github.com/stefpeschel/NetCoMi; see FANG et al. (2017) <doi:10.1089/cmb.2017.0054> for the algorithm details
- 'FlashWeave'
- FlashWeavenetwork; Local-to-global learning framework; belong to the category of conditional dependence and graphical models; good performance on heterogenous datasets to find direct associations among taxa; see https://github.com/meringlab/FlashWeave.jl for installing- julialanguage and- FlashWeavepackage; julia must be in the computer system env path, otherwise the program can not find it; see Tackmann et al. (2019) <doi:10.1016/j.cels.2019.08.002> for the algorithm details
- 'beemStatic'
- beemStaticnetwork; extend generalized Lotka-Volterra model to cases of cross-sectional datasets to infer interaction among taxa based on expectation-maximization algorithm; see https://github.com/CSB5/BEEM-static for installing the R package; see Li et al. (2021) <doi:10.1371/journal.pcbi.1009343> for the algorithm details
 
- COR_p_thres
- default 0.01; the p value threshold for the correlation-based network. 
- COR_p_adjust
- default "fdr"; p value adjustment method, see - methodparameter of- p.adjustfunction for available options, in which- COR_p_adjust = "none"means giving up the p value adjustment.
- COR_weight
- default TRUE; whether use correlation coefficient as the weight of edges; FALSE represents weight = 1 for all edges. 
- COR_cut
- default 0.6; correlation coefficient threshold for the correlation network. 
- COR_optimization
- default FALSE; whether use random matrix theory (RMT) based method to determine the correlation coefficient; see https://doi.org/10.1186/1471-2105-13-113 
- COR_optimization_low_high
- default - c(0.01, 0.8); the low and high value threshold used for the RMT optimization; only useful when COR_optimization = TRUE.
- COR_optimization_seq
- default 0.01; the interval of correlation coefficient used for RMT optimization; only useful when COR_optimization = TRUE. 
- SpiecEasi_method
- default "mb"; either 'glasso' or 'mb';see spiec.easi function in package SpiecEasi and https://github.com/zdk123/SpiecEasi. 
- FlashWeave_tempdir
- default NULL; The temporary directory used to save the temporary files for running FlashWeave; If not assigned, use the system user temp. 
- FlashWeave_meta_data
- default FALSE; whether use env data for the optimization, If TRUE, the function automatically find the - env_datain the object and generate a file for meta_data_path parameter of FlashWeave package.
- FlashWeave_other_para
- default - "alpha=0.01,sensitive=true,heterogeneous=true"; the parameters passed to julia FlashWeave package; user can change the parameters or add more according to FlashWeave help document; An exception is meta_data_path parameter as it is generated based on the data inside the object, see FlashWeave_meta_data parameter for the description.
- beemStatic_t_strength
- default 0.001; for network_method = "beemStatic"; the threshold used to limit the number of interactions (strength); same with the t.strength parameter in showInteraction function of beemStatic package. 
- beemStatic_t_stab
- default 0.8; for network_method = "beemStatic"; the threshold used to limit the number of interactions (stability); same with the t.stab parameter in showInteraction function of beemStatic package. 
- add_taxa_name
- default "Phylum"; one or more taxonomic rank name; used to add taxonomic rank name to network node properties. 
- delete_unlinked_nodes
- default TRUE; whether delete the nodes without any link. 
- usename_rawtaxa_when_taxalevel_notOTU
- default FALSE; whether use OTU name as representatives of taxa when - taxa_level != "OTU". Default- FALSEmeans using taxonomic information of- taxa_levelinstead of OTU name.
- ...
- parameters pass to - SpiecEasi::spiec.easiwhen- network_method = "SpiecEasi"; pass to- NetCoMi::netConstructwhen- network_method = "gcoda"; pass to- beemStatic::func.EMwhen- network_method = "beemStatic".
Returns
res_network stored in object.
Examples
\dontrun{
# for correlation network
t1 <- trans_network$new(dataset = dataset, cor_method = "pearson", 
		taxa_level = "OTU", filter_thres = 0.001)
t1$cal_network(COR_p_thres = 0.05, COR_cut = 0.6)
t1 <- trans_network$new(dataset = dataset, cor_method = NULL, filter_thres = 0.003)
t1$cal_network(network_method = "SpiecEasi", SpiecEasi_method = "mb")
t1 <- trans_network$new(dataset = dataset, cor_method = NULL, filter_thres = 0.005)
t1$cal_network(network_method = "beemStatic")
t1 <- trans_network$new(dataset = dataset, cor_method = NULL, filter_thres = 0.001)
t1$cal_network(network_method = "FlashWeave")
}
Method cal_module()
Calculate network modules and add module names to the network node properties.
Usage
trans_network$cal_module( method = "cluster_fast_greedy", module_name_prefix = "M" )
Arguments
- method
- default "cluster_fast_greedy"; the method used to find the optimal community structure of a graph; the following are available functions (options) from igraph package: 
 - "cluster_fast_greedy",- "cluster_walktrap",- "cluster_edge_betweenness",
 - "cluster_infomap",- "cluster_label_prop",- "cluster_leading_eigen",
 - "cluster_louvain",- "cluster_spinglass",- "cluster_optimal".
 For the details of these functions, please see the help document, such as- help(cluster_fast_greedy); Note that the default- "cluster_fast_greedy"method can not be applied to directed network. If directed network is provided, the function can automatically switch the default method from- "cluster_fast_greedy"to- "cluster_walktrap".
- module_name_prefix
- default "M"; the prefix of module names; module names are made of the module_name_prefix and numbers; numbers are assigned according to the sorting result of node numbers in modules with decreasing trend. 
Returns
res_network with modules, stored in object.
Examples
\donttest{
t1 <- trans_network$new(dataset = dataset, cor_method = "pearson", 
		taxa_level = "OTU", filter_thres = 0.0002)
t1$cal_network(COR_p_thres = 0.01, COR_cut = 0.6)
t1$cal_module(method = "cluster_fast_greedy")
}
Method save_network()
Save network as gexf style, which can be opened by Gephi (https://gephi.org/).
Usage
trans_network$save_network(filepath = "network.gexf")
Arguments
- filepath
- default "network.gexf"; file path to save the network. 
Returns
None
Examples
\dontrun{
t1$save_network(filepath = "network.gexf")
}
Method cal_network_attr()
Calculate network properties.
Usage
trans_network$cal_network_attr()
Returns
res_network_attr stored in object.
Examples
\donttest{
t1$cal_network_attr()
}
Method get_node_table()
Get the node property table. The properties include the node names, modules allocation, degree, betweenness, abundance, taxonomy, within-module connectivity (zi) and among-module connectivity (Pi) <doi:10.1186/1471-2105-13-113; 10.1016/j.geoderma.2022.115866>.
Usage
trans_network$get_node_table(node_roles = TRUE)
Arguments
- node_roles
- default TRUE; whether calculate the node roles <doi:10.1038/nature03288; 10.1186/1471-2105-13-113>. The role of node i is characterized by its within-module connectivity (zi) and among-module connectivity (Pi) as follows - z_i = \dfrac{k_{ib} - \bar{k_b}}{\sigma_{k_b}}- P_i = 1 - \displaystyle\sum_{c=1}^{N_M} \biggl(\frac{k_{ic}}{k_i}\biggr)^2- where - k_{ib}is the number of links of node- ito other nodes in its module- b,- \bar{k_b}and- \sigma_{k_b}are the average and standard deviation of within-module connectivity, respectively over all the nodes in module- b,- k_iis the number of links of node- iin the whole network,- k_{ic}is the number of links from node- ito nodes in module- c, and- N_Mis the number of modules in the network.
Returns
res_node_table in object; Abundance expressed as a percentage; 
betweenness_centrality: betweenness centrality; betweenness_centrality: closeness centrality; eigenvector_centrality: eigenvector centrality; 
z: within-module connectivity; p: among-module connectivity.
Examples
\donttest{
t1$get_node_table(node_roles = TRUE)
}
Method get_edge_table()
Get the edge property table, including connected nodes, label and weight.
Usage
trans_network$get_edge_table()
Returns
res_edge_table in object.
Examples
\donttest{
t1$get_edge_table()
}
Method get_adjacency_matrix()
Get the adjacency matrix from the network graph.
Usage
trans_network$get_adjacency_matrix(...)
Arguments
- ...
- parameters passed to as_adjacency_matrix function of - igraphpackage.
Returns
res_adjacency_matrix in object.
Examples
\donttest{
t1$get_adjacency_matrix(attr = "weight")
}
Method plot_network()
Plot the network based on a series of methods from other packages, such as igraph, ggraph and networkD3. 
The networkD3 package provides dynamic network. It is especially useful for a glimpse of the whole network structure and finding 
the interested nodes and edges in a large network. In contrast, the igraph and ggraph methods are suitable for relatively small network.
Usage
trans_network$plot_network(
  method = c("igraph", "ggraph", "networkD3")[1],
  node_label = "name",
  node_color = NULL,
  ggraph_layout = "fr",
  ggraph_node_size = 2,
  ggraph_node_text = TRUE,
  ggraph_text_color = NULL,
  ggraph_text_size = 3,
  networkD3_node_legend = TRUE,
  networkD3_zoom = TRUE,
  ...
)Arguments
- method
- default "igraph"; The available options: - 'igraph'
- call - plot.igraphfunction in- igraphpackage for a static network; see plot.igraph for the parameters
- 'ggraph'
- call - ggraphfunction in- ggraphpackage for a static network
- 'networkD3'
- use forceNetwork function in - networkD3package for a dynamic network; see forceNetwork function for the parameters
 
- node_label
- default "name"; node label shown in the plot for - method = "ggraph"or- method = "networkD3"; Please see the column names of object$res_node_table, which is the returned table of function object$get_node_table; User can select other column names in res_node_table.
- node_color
- default NULL; node color assignment for - method = "ggraph"or- method = "networkD3"; Select a column name of- object$res_node_table, such as "module".
- ggraph_layout
- default "fr"; for - method = "ggraph"; see- layoutparameter of- create_layoutfunction in- ggraphpackage.
- ggraph_node_size
- default 2; for - method = "ggraph"; the node size.
- ggraph_node_text
- default TRUE; for - method = "ggraph"; whether show the label text of nodes.
- ggraph_text_color
- default NULL; for - method = "ggraph"; a column name of object$res_node_table used to assign label text colors.
- ggraph_text_size
- default 3; for - method = "ggraph"; the node label text size.
- networkD3_node_legend
- default TRUE; used for - method = "networkD3"; logical value to enable node colour legends; Please see the legend parameter in networkD3::forceNetwork function.
- networkD3_zoom
- default TRUE; used for - method = "networkD3"; logical value to enable (TRUE) or disable (FALSE) zooming; Please see the zoom parameter in networkD3::forceNetwork function.
- ...
- parameters passed to - plot.igraphfunction when- method = "igraph"or forceNetwork function when- method = "networkD3".
Returns
network plot.
Examples
\donttest{
t1$plot_network(method = "igraph", layout = layout_with_kk)
t1$plot_network(method = "ggraph", node_color = "module")
t1$plot_network(method = "networkD3", node_color = "module")
}
Method cal_eigen()
Calculate eigengenes of modules, i.e. the first principal component based on PCA analysis, and the percentage of variance <doi:10.1186/1471-2105-13-113>.
Usage
trans_network$cal_eigen()
Returns
res_eigen and res_eigen_expla in object.
Examples
\donttest{
t1$cal_eigen()
}
Method plot_taxa_roles()
Plot the roles or metrics of nodes based on the res_node_table data (coming from function get_node_table) stored in the object.
Usage
trans_network$plot_taxa_roles(
  use_type = c(1, 2)[1],
  roles_color_background = FALSE,
  roles_color_values = NULL,
  add_label = FALSE,
  add_label_group = "Network hubs",
  add_label_text = "name",
  label_text_size = 4,
  label_text_color = "grey50",
  label_text_italic = FALSE,
  label_text_parse = FALSE,
  plot_module = FALSE,
  x_lim = c(0, 1),
  use_level = "Phylum",
  show_value = c("z", "p"),
  show_number = 1:10,
  plot_color = "Phylum",
  plot_shape = "taxa_roles",
  plot_size = "Abundance",
  color_values = RColorBrewer::brewer.pal(12, "Paired"),
  shape_values = c(16, 17, 7, 8, 15, 18, 11, 10, 12, 13, 9, 3, 4, 0, 1, 2, 14),
  ...
)Arguments
- use_type
- default 1; 1 or 2; 1 represents taxa roles plot (node roles include Module hubs, Network hubs, Connectors and Peripherals <doi:10.1038/nature03288; 10.1186/1471-2105-13-113>); 2 represents the layered plot with taxa as x axis and the index (e.g., Zi and Pi) as y axis. Please refer to - res_node_tabledata stored in the object for the detailed information.
- roles_color_background
- default FALSE; for use_type=1; TRUE: use background colors for each area; FALSE: use classic point colors. 
- roles_color_values
- default NULL; for use_type=1; color palette for background or points. 
- add_label
- default FALSE; for use_type = 1; whether add labels for the points. 
- add_label_group
- default "Network hubs"; If add_label = TRUE; which part of tax_roles is used to show labels; character vectors. 
- add_label_text
- default "name"; If add_label = TRUE; which column of object$res_node_table is used to label the text. 
- label_text_size
- default 4; The text size of the label. 
- label_text_color
- default "grey50"; The text color of the label. 
- label_text_italic
- default FALSE; whether use italic style for the label text. 
- label_text_parse
- default FALSE; whether parse the label text. See the parse parameter in - ggrepel::geom_text_repelfunction.
- plot_module
- default FALSE; for use_type=1; whether plot the modules information. 
- x_lim
- default c(0, 1); for use_type=1; x axis range when roles_color_background = FALSE. 
- use_level
- default "Phylum"; for use_type=2; used taxonomic level in x axis. 
- show_value
- default c("z", "p"); for use_type=2; indexes used in y axis. Please see - res_node_tablein the object for other available indexes.
- show_number
- default 1:10; for use_type=2; showed number in x axis, sorting according to the nodes number. 
- plot_color
- default "Phylum"; for use_type=2; variable for color. 
- plot_shape
- default "taxa_roles"; for use_type=2; variable for shape. 
- plot_size
- default "Abundance"; for use_type=2; used for point size; a fixed number (e.g. 5) is also acceptable. 
- color_values
- default RColorBrewer::brewer.pal(12, "Paired"); for use_type=2; color vector. 
- shape_values
- default c(16, 17, 7, 8, 15, 18, 11, 10, 12, 13, 9, 3, 4, 0, 1, 2, 14); for use_type=2; shape vector, see ggplot2 tutorial for the shape meaning. 
- ...
- parameters pass to - geom_pointfunction of ggplot2 package.
Returns
ggplot.
Examples
\donttest{
t1$plot_taxa_roles(roles_color_background = FALSE)
}
Method subset_network()
Subset of the network.
Usage
trans_network$subset_network(node = NULL, edge = NULL, rm_single = TRUE)
Arguments
- node
- default NULL; provide the node names that you want to use in the sub-network. 
- edge
- default NULL; provide the edge name needed; must be one of "+" or "-". 
- rm_single
- default TRUE; whether remove the nodes without any edge in the sub-network. So this function can also be used to remove the nodes withou any edge when node and edge are both NULL. 
Returns
a new network
Examples
\donttest{
t1$subset_network(node = t1$res_node_table %>% base::subset(module == "M1") %>% 
  rownames, rm_single = TRUE)
# return a sub network that contains all nodes of module M1
}
Method cal_powerlaw()
Fit degrees to a power law distribution. First, perform a bootstrapping hypothesis test to determine whether degrees follow a power law distribution. If the distribution follows power law, then fit degrees to power law distribution and return the parameters.
Usage
trans_network$cal_powerlaw(...)
Arguments
- ...
- parameters pass to bootstrap_p function in poweRlaw package. 
Returns
res_powerlaw_p and res_powerlaw_fit; see poweRlaw::bootstrap_p function for the bootstrapping p value details;
see igraph::fit_power_law function for the power law fit return details.
Examples
\donttest{
t1$cal_powerlaw()
}
Method cal_sum_links()
This function is used to sum the links number from one taxa to another or in the same taxa, for example, at Phylum level. This is very useful to fast see how many nodes are connected between different taxa or within the taxa.
Usage
trans_network$cal_sum_links(taxa_level = "Phylum")
Arguments
- taxa_level
- default "Phylum"; taxonomic rank. 
Returns
res_sum_links_pos and res_sum_links_neg in object.
Examples
\donttest{
t1$cal_sum_links(taxa_level = "Phylum")
}
Method plot_sum_links()
Plot the summed linkages among taxa.
Usage
trans_network$plot_sum_links(
  plot_pos = TRUE,
  plot_num = NULL,
  color_values = RColorBrewer::brewer.pal(8, "Dark2"),
  method = c("chorddiag", "circlize")[1],
  ...
)Arguments
- plot_pos
- default TRUE; If TRUE, plot the summed positive linkages; If FALSE, plot the summed negative linkages. 
- plot_num
- default NULL; number of taxa presented in the plot. 
- color_values
- default RColorBrewer::brewer.pal(8, "Dark2"); colors palette for taxa. 
- method
- default c("chorddiag", "circlize")[1]; chorddiag package <https://github.com/mattflor/chorddiag> or circlize package. 
- ...
- pass to - chorddiag::chorddiagfunction when- method = "chorddiag"or- circlize::chordDiagramfunction when- method = "circlize". Note that for- circlize::chordDiagramfunction,- keep.diagonal,- symmetricand- self.linkparameters have been fixed to fit the input data.
Returns
please see the invoked function.
Examples
\dontrun{
test1$plot_sum_links(method = "chorddiag", plot_pos = TRUE, plot_num = 10)
test1$plot_sum_links(method = "circlize", transparency = 0.2, 
	  annotationTrackHeight = circlize::mm_h(c(5, 5)))
}
Method random_network()
Generate random networks, compare them with the empirical network and get the p value of topological properties.
The generation of random graph is based on the erdos.renyi.game function of igraph package.
The numbers of vertices and edges in the random graph are same with the empirical network stored in the object.
Usage
trans_network$random_network(runs = 100, output_sim = FALSE)
Arguments
- runs
- default 100; simulation number of random network. 
- output_sim
- default FALSE; whether output each simulated network result. 
Returns
a data.frame with the following components:
- Observed
- Topological properties of empirical network 
- Mean_sim
- Mean of properties of simulated networks 
- SD_sim
- SD of properties of simulated networks 
- p_value
- Significance, i.e. p values 
When output_sim = TRUE, the columns from the five to the last are each simulated result.
Examples
\dontrun{
t1$random_network(runs = 100)
}
Method trans_comm()
Transform classifed features to community-like microtable object for further analysis, such as module-taxa table.
Usage
trans_network$trans_comm(use_col = "module", abundance = TRUE)
Arguments
- use_col
- default "module"; which column to use as the 'community'; must be one of the name of res_node_table from function - get_node_table.
- abundance
- default TRUE; whether sum abundance of taxa. TRUE: sum the abundance for a taxon across all samples; FALSE: sum the frequency for a taxon across all samples. 
Returns
a new microtable class.
Examples
\donttest{
t2 <- t1$trans_comm(use_col = "module")
}
Method print()
Print the trans_network object.
Usage
trans_network$print()
Method clone()
The objects of this class are cloneable with this method.
Usage
trans_network$clone(deep = FALSE)
Arguments
- deep
- Whether to make a deep clone. 
Examples
## ------------------------------------------------
## Method `trans_network$new`
## ------------------------------------------------
data(dataset)
# for correlation network
t1 <- trans_network$new(dataset = dataset, cor_method = "pearson", 
		taxa_level = "OTU", filter_thres = 0.0002)
# for non-correlation network
t1 <- trans_network$new(dataset = dataset, cor_method = NULL)
## ------------------------------------------------
## Method `trans_network$cal_network`
## ------------------------------------------------
## Not run: 
# for correlation network
t1 <- trans_network$new(dataset = dataset, cor_method = "pearson", 
		taxa_level = "OTU", filter_thres = 0.001)
t1$cal_network(COR_p_thres = 0.05, COR_cut = 0.6)
t1 <- trans_network$new(dataset = dataset, cor_method = NULL, filter_thres = 0.003)
t1$cal_network(network_method = "SpiecEasi", SpiecEasi_method = "mb")
t1 <- trans_network$new(dataset = dataset, cor_method = NULL, filter_thres = 0.005)
t1$cal_network(network_method = "beemStatic")
t1 <- trans_network$new(dataset = dataset, cor_method = NULL, filter_thres = 0.001)
t1$cal_network(network_method = "FlashWeave")
## End(Not run)
## ------------------------------------------------
## Method `trans_network$cal_module`
## ------------------------------------------------
t1 <- trans_network$new(dataset = dataset, cor_method = "pearson", 
		taxa_level = "OTU", filter_thres = 0.0002)
t1$cal_network(COR_p_thres = 0.01, COR_cut = 0.6)
t1$cal_module(method = "cluster_fast_greedy")
## ------------------------------------------------
## Method `trans_network$save_network`
## ------------------------------------------------
## Not run: 
t1$save_network(filepath = "network.gexf")
## End(Not run)
## ------------------------------------------------
## Method `trans_network$cal_network_attr`
## ------------------------------------------------
t1$cal_network_attr()
## ------------------------------------------------
## Method `trans_network$get_node_table`
## ------------------------------------------------
t1$get_node_table(node_roles = TRUE)
## ------------------------------------------------
## Method `trans_network$get_edge_table`
## ------------------------------------------------
t1$get_edge_table()
## ------------------------------------------------
## Method `trans_network$get_adjacency_matrix`
## ------------------------------------------------
t1$get_adjacency_matrix(attr = "weight")
## ------------------------------------------------
## Method `trans_network$plot_network`
## ------------------------------------------------
t1$plot_network(method = "igraph", layout = layout_with_kk)
t1$plot_network(method = "ggraph", node_color = "module")
t1$plot_network(method = "networkD3", node_color = "module")
## ------------------------------------------------
## Method `trans_network$cal_eigen`
## ------------------------------------------------
t1$cal_eigen()
## ------------------------------------------------
## Method `trans_network$plot_taxa_roles`
## ------------------------------------------------
t1$plot_taxa_roles(roles_color_background = FALSE)
## ------------------------------------------------
## Method `trans_network$subset_network`
## ------------------------------------------------
t1$subset_network(node = t1$res_node_table %>% base::subset(module == "M1") %>% 
  rownames, rm_single = TRUE)
# return a sub network that contains all nodes of module M1
## ------------------------------------------------
## Method `trans_network$cal_powerlaw`
## ------------------------------------------------
t1$cal_powerlaw()
## ------------------------------------------------
## Method `trans_network$cal_sum_links`
## ------------------------------------------------
t1$cal_sum_links(taxa_level = "Phylum")
## ------------------------------------------------
## Method `trans_network$plot_sum_links`
## ------------------------------------------------
## Not run: 
test1$plot_sum_links(method = "chorddiag", plot_pos = TRUE, plot_num = 10)
test1$plot_sum_links(method = "circlize", transparency = 0.2, 
	  annotationTrackHeight = circlize::mm_h(c(5, 5)))
## End(Not run)
## ------------------------------------------------
## Method `trans_network$random_network`
## ------------------------------------------------
## Not run: 
t1$random_network(runs = 100)
## End(Not run)
## ------------------------------------------------
## Method `trans_network$trans_comm`
## ------------------------------------------------
t2 <- t1$trans_comm(use_col = "module")