| Title: | Inferring Cell-Specific Gene Regulatory Network |
|---|---|
| Description: | Inferring cell-type specific gene regulatory network using sparse regression model from single-cell transcriptomic data. |
| Authors: | Meng Xu [aut, cre] (ORCID: <https://orcid.org/0000-0002-8300-1054>) |
| Maintainer: | Meng Xu <[email protected]> |
| License: | MIT + file LICENSE |
| Version: | 1.2.3 |
| Built: | 2026-05-24 15:12:22 UTC |
| Source: | https://github.com/mengxu98/infercsn |
Inferring cell-type specific gene regulatory network using sparse regression model from single-cell transcriptomic data.
Meng Xu (Maintainer), [email protected]
https://mengxu98.github.io/inferCSN/
Useful links:
Calculates accuracy metric
calculate_accuracy(network_table, ground_truth)calculate_accuracy(network_table, ground_truth)
network_table |
A data frame of predicted network structure |
ground_truth |
A data frame of ground truth network |
A list containing the metric
data(example_matrix) data("example_ground_truth") network_table <- inferCSN(example_matrix) calculate_accuracy( network_table, example_ground_truth )data(example_matrix) data("example_ground_truth") network_table <- inferCSN(example_matrix) calculate_accuracy( network_table, example_ground_truth )
Calculates AUROC and AUPRC metrics with optional visualization
calculate_auc( network_table, ground_truth, return_plot = FALSE, line_color = "#1563cc", line_width = 1, tag_levels = "A" )calculate_auc( network_table, ground_truth, return_plot = FALSE, line_color = "#1563cc", line_width = 1, tag_levels = "A" )
network_table |
A data frame of predicted network structure |
ground_truth |
A data frame of ground truth network |
return_plot |
Logical value indicating whether to generate plots |
line_color |
Color for plot lines |
line_width |
Width for plot lines |
tag_levels |
Tag levels for plot annotations |
A list containing metrics and optional plots
data(example_matrix) data("example_ground_truth") network_table <- inferCSN(example_matrix) calculate_auc( network_table, example_ground_truth, return_plot = TRUE )data(example_matrix) data("example_ground_truth") network_table <- inferCSN(example_matrix) calculate_auc( network_table, example_ground_truth, return_plot = TRUE )
Calculates AUPRC metric with optional visualization
calculate_auprc( network_table, ground_truth, return_plot = FALSE, line_color = "#1563cc", line_width = 1 )calculate_auprc( network_table, ground_truth, return_plot = FALSE, line_color = "#1563cc", line_width = 1 )
network_table |
A data frame of predicted network structure |
ground_truth |
A data frame of ground truth network |
return_plot |
Logical value indicating whether to generate plot |
line_color |
Color for plot lines |
line_width |
Width for plot lines |
A list containing metric and optional plot
data(example_matrix) data("example_ground_truth") network_table <- inferCSN(example_matrix) calculate_auprc( network_table, example_ground_truth, return_plot = TRUE )data(example_matrix) data("example_ground_truth") network_table <- inferCSN(example_matrix) calculate_auprc( network_table, example_ground_truth, return_plot = TRUE )
Calculates AUROC metric with optional visualization
calculate_auroc( network_table, ground_truth, return_plot = FALSE, line_color = "#1563cc", line_width = 1 )calculate_auroc( network_table, ground_truth, return_plot = FALSE, line_color = "#1563cc", line_width = 1 )
network_table |
A data frame of predicted network structure |
ground_truth |
A data frame of ground truth network |
return_plot |
Logical value indicating whether to generate plot |
line_color |
Color for plot lines |
line_width |
Width for plot lines |
A list containing metric and optional plot
data(example_matrix) data("example_ground_truth") network_table <- inferCSN(example_matrix) calculate_auroc( network_table, example_ground_truth, return_plot = TRUE )data(example_matrix) data("example_ground_truth") network_table <- inferCSN(example_matrix) calculate_auroc( network_table, example_ground_truth, return_plot = TRUE )
Calculates F1 score
calculate_f1(network_table, ground_truth)calculate_f1(network_table, ground_truth)
network_table |
A data frame of predicted network structure |
ground_truth |
A data frame of ground truth network |
A list containing the metric
data(example_matrix) data("example_ground_truth") network_table <- inferCSN(example_matrix) calculate_f1( network_table, example_ground_truth )data(example_matrix) data("example_ground_truth") network_table <- inferCSN(example_matrix) calculate_f1( network_table, example_ground_truth )
Rank TFs and genes in network
calculate_gene_rank( network_table, regulators = NULL, targets = NULL, directed = FALSE )calculate_gene_rank( network_table, regulators = NULL, targets = NULL, directed = FALSE )
network_table |
The weight data table of network. |
regulators |
Regulators list. |
targets |
Targets list. |
directed |
Whether the network is directed. |
A table of gene rank.
data(example_matrix) network_table <- inferCSN(example_matrix) head(calculate_gene_rank(network_table)) head(calculate_gene_rank(network_table, regulators = "g1")) head(calculate_gene_rank(network_table, targets = "g1"))data(example_matrix) network_table <- inferCSN(example_matrix) head(calculate_gene_rank(network_table)) head(calculate_gene_rank(network_table, regulators = "g1")) head(calculate_gene_rank(network_table, targets = "g1"))
Calculates Jaccard Index (JI) metric
calculate_ji(network_table, ground_truth)calculate_ji(network_table, ground_truth)
network_table |
A data frame of predicted network structure |
ground_truth |
A data frame of ground truth network |
A list containing the metric
data(example_matrix) data("example_ground_truth") network_table <- inferCSN(example_matrix) calculate_ji( network_table, example_ground_truth )data(example_matrix) data("example_ground_truth") network_table <- inferCSN(example_matrix) calculate_ji( network_table, example_ground_truth )
Calculates comprehensive performance metrics for evaluating predicted network structures, including classification performance, precision-recall metrics, and network topology metrics.
calculate_metrics( network_table, ground_truth, metric_type = c("all", "auc", "auroc", "auprc", "precision", "recall", "f1", "accuracy", "si", "ji"), return_plot = FALSE, line_color = "#1563cc", line_width = 1 )calculate_metrics( network_table, ground_truth, metric_type = c("all", "auc", "auroc", "auprc", "precision", "recall", "f1", "accuracy", "si", "ji"), return_plot = FALSE, line_color = "#1563cc", line_width = 1 )
network_table |
A data frame. |
ground_truth |
A data frame of ground truth network with the same format as network_table. |
metric_type |
The type of metric to return.
Default is
|
return_plot |
Whether to generate visualization plots.
Default is |
line_color |
Color for plot lines.
Default is |
line_width |
Width for plot lines.
Default is |
A list containing:
metrics A data frame with requested metrics.
plot A plot object if return_plot = TRUE (optional).
data(example_matrix) data("example_ground_truth") network_table <- inferCSN(example_matrix) calculate_metrics( network_table, example_ground_truth, return_plot = TRUE ) calculate_metrics( network_table, example_ground_truth, metric_type = "auroc" )data(example_matrix) data("example_ground_truth") network_table <- inferCSN(example_matrix) calculate_metrics( network_table, example_ground_truth, return_plot = TRUE ) calculate_metrics( network_table, example_ground_truth, metric_type = "auroc" )
Calculates precision metric
calculate_precision(network_table, ground_truth)calculate_precision(network_table, ground_truth)
network_table |
A data frame of predicted network structure |
ground_truth |
A data frame of ground truth network |
A list containing the metric
data(example_matrix) data("example_ground_truth") network_table <- inferCSN(example_matrix) calculate_precision( network_table, example_ground_truth )data(example_matrix) data("example_ground_truth") network_table <- inferCSN(example_matrix) calculate_precision( network_table, example_ground_truth )
Calculates recall metric
calculate_recall(network_table, ground_truth)calculate_recall(network_table, ground_truth)
network_table |
A data frame of predicted network structure |
ground_truth |
A data frame of ground truth network |
A list containing the metric
data(example_matrix) data("example_ground_truth") network_table <- inferCSN(example_matrix) calculate_recall( network_table, example_ground_truth )data(example_matrix) data("example_ground_truth") network_table <- inferCSN(example_matrix) calculate_recall( network_table, example_ground_truth )
Calculates Set Intersection (SI) metric
calculate_si(network_table, ground_truth)calculate_si(network_table, ground_truth)
network_table |
A data frame of predicted network structure |
ground_truth |
A data frame of ground truth network |
A list containing the metric
data(example_matrix) data("example_ground_truth") network_table <- inferCSN(example_matrix) calculate_si( network_table, example_ground_truth )data(example_matrix) data("example_ground_truth") network_table <- inferCSN(example_matrix) calculate_si( network_table, example_ground_truth )
Filter and sort matrix
filter_sort_matrix(network_matrix, regulators = NULL, targets = NULL)filter_sort_matrix(network_matrix, regulators = NULL, targets = NULL)
network_matrix |
The matrix of network weight. |
regulators |
Regulators list. |
targets |
Targets list. |
Filtered and sorted matrix
data(example_matrix) network_table <- inferCSN(example_matrix) colnames(network_table) <- c("row", "col", "value") network_matrix <- thisutils::table_to_matrix(network_table) filter_sort_matrix(network_matrix)[1:6, 1:6] filter_sort_matrix( network_matrix, regulators = c("g1", "g2"), targets = c("g3", "g4") )data(example_matrix) network_table <- inferCSN(example_matrix) colnames(network_table) <- c("row", "col", "value") network_matrix <- thisutils::table_to_matrix(network_table) filter_sort_matrix(network_matrix)[1:6, 1:6] filter_sort_matrix( network_matrix, regulators = c("g1", "g2"), targets = c("g3", "g4") )
Sparse regression model
fit_srm( x, y, cross_validation = FALSE, seed = 1, penalty = "L0", regulators_num = ncol(x), n_folds = 5, verbose = TRUE, ... )fit_srm( x, y, cross_validation = FALSE, seed = 1, penalty = "L0", regulators_num = ncol(x), n_folds = 5, verbose = TRUE, ... )
x |
The matrix of regulators. |
y |
The vector of target. |
cross_validation |
Whether to use cross-validation.
Default is |
seed |
The random seed for cross-validation.
Default is |
penalty |
The type of regularization.
This can take either one of the following choices: |
regulators_num |
The number of regulators for target. |
n_folds |
The number of folds for cross-validation.
Default is |
verbose |
Whether to print progress messages.
Default is |
... |
Parameters for other methods. |
A list of the sparse regression model. The list has the three components: model, metrics, and coefficients.
data(example_matrix) fit_srm( x = example_matrix[, -1], y = example_matrix[, 1] )data(example_matrix) fit_srm( x = example_matrix[, -1], y = example_matrix[, 1] )
inferring cell-type specific gene regulatory network
inferCSN( object, penalty = c("L0", "L0L1", "L0L2"), cross_validation = FALSE, seed = 1, n_folds = 5, subsampling_method = c("sample", "meta_cells", "pseudobulk"), subsampling_ratio = 1, r_squared_threshold = 0, regulators = NULL, targets = NULL, cores = 1, verbose = TRUE, ... ) ## S4 method for signature 'matrix' inferCSN( object, penalty = c("L0", "L0L1", "L0L2"), cross_validation = FALSE, seed = 1, n_folds = 5, subsampling_method = c("sample", "meta_cells", "pseudobulk"), subsampling_ratio = 1, r_squared_threshold = 0, regulators = NULL, targets = NULL, cores = 1, verbose = TRUE, ... ) ## S4 method for signature 'sparseMatrix' inferCSN( object, penalty = c("L0", "L0L1", "L0L2"), cross_validation = FALSE, seed = 1, n_folds = 5, subsampling_method = c("sample", "meta_cells", "pseudobulk"), subsampling_ratio = 1, r_squared_threshold = 0, regulators = NULL, targets = NULL, cores = 1, verbose = TRUE, ... )inferCSN( object, penalty = c("L0", "L0L1", "L0L2"), cross_validation = FALSE, seed = 1, n_folds = 5, subsampling_method = c("sample", "meta_cells", "pseudobulk"), subsampling_ratio = 1, r_squared_threshold = 0, regulators = NULL, targets = NULL, cores = 1, verbose = TRUE, ... ) ## S4 method for signature 'matrix' inferCSN( object, penalty = c("L0", "L0L1", "L0L2"), cross_validation = FALSE, seed = 1, n_folds = 5, subsampling_method = c("sample", "meta_cells", "pseudobulk"), subsampling_ratio = 1, r_squared_threshold = 0, regulators = NULL, targets = NULL, cores = 1, verbose = TRUE, ... ) ## S4 method for signature 'sparseMatrix' inferCSN( object, penalty = c("L0", "L0L1", "L0L2"), cross_validation = FALSE, seed = 1, n_folds = 5, subsampling_method = c("sample", "meta_cells", "pseudobulk"), subsampling_ratio = 1, r_squared_threshold = 0, regulators = NULL, targets = NULL, cores = 1, verbose = TRUE, ... )
object |
The input data for inferring network. |
penalty |
The type of regularization.
This can take either one of the following choices: |
cross_validation |
Whether to use cross-validation.
Default is |
seed |
The random seed for cross-validation.
Default is |
n_folds |
The number of folds for cross-validation.
Default is |
subsampling_method |
The method to use for subsampling.
Options are |
subsampling_ratio |
The percent of all samples used for fit_srm.
Default is |
r_squared_threshold |
Threshold of R² coefficient.
Default is |
regulators |
The regulator genes for which to infer the regulatory network. |
targets |
The target genes for which to infer the regulatory network. |
cores |
The number of cores to use for parallelization with foreach::foreach.
Default is |
verbose |
Whether to print progress messages.
Default is |
... |
Parameters for other methods. |
A data table of regulator-target regulatory relationships,
which has three columns: regulator, target, and weight.
data(example_matrix) network_table <- inferCSN( example_matrix ) head(network_table) inferCSN( example_matrix, regulators = c("g1", "g2"), targets = c("g3", "g4") ) inferCSN( example_matrix, regulators = c("g1", "g2"), targets = c("g3", "g0") )data(example_matrix) network_table <- inferCSN( example_matrix ) head(network_table) inferCSN( example_matrix, regulators = c("g1", "g2"), targets = c("g3", "g4") ) inferCSN( example_matrix, regulators = c("g1", "g2"), targets = c("g3", "g0") )
The inferCSN logo, using ASCII or Unicode characters Use cli::ansi_strip to get rid of the colors.
infercsn_logo(unicode = cli::is_utf8_output())infercsn_logo(unicode = cli::is_utf8_output())
unicode |
Unicode symbols on UTF-8 platforms. Default is cli::is_utf8_output. |
https://github.com/tidyverse/tidyverse/blob/main/R/logo.R
infercsn_logo()infercsn_logo()
This function detects metacells from a single-cell gene expression matrix using dimensionality reduction and clustering techniques.
meta_cells( matrix, genes_use = NULL, genes_exclude = NULL, var_genes_num = min(1000, nrow(matrix)), gamma = 10, knn_k = 5, do_scale = TRUE, pc_num = 25, fast_pca = FALSE, do_approx = FALSE, approx_num = 20000, directed = FALSE, use_nn2 = TRUE, seed = 1, cluster_method = "walktrap", block_size = 10000, weights = NULL, do_median_norm = FALSE, ... )meta_cells( matrix, genes_use = NULL, genes_exclude = NULL, var_genes_num = min(1000, nrow(matrix)), gamma = 10, knn_k = 5, do_scale = TRUE, pc_num = 25, fast_pca = FALSE, do_approx = FALSE, approx_num = 20000, directed = FALSE, use_nn2 = TRUE, seed = 1, cluster_method = "walktrap", block_size = 10000, weights = NULL, do_median_norm = FALSE, ... )
matrix |
A gene expression matrix where rows represent genes and columns represent cells. |
genes_use |
A character vector specifying genes to be used for PCA dimensionality reduction.
Default is |
genes_exclude |
A character vector specifying genes to be excluded from PCA computation.
Default is |
var_genes_num |
Number of most variable genes to select when |
gamma |
Default is |
knn_k |
Default is |
do_scale |
Whether to standardize (center and scale) gene expression values before PCA.
Default is |
pc_num |
Default is |
fast_pca |
Default is |
do_approx |
Default is |
approx_num |
Default is |
directed |
Default is |
use_nn2 |
Default is |
seed |
Default is |
cluster_method |
Default is |
block_size |
Default is |
weights |
Default is |
do_median_norm |
Default is |
... |
Additional arguments passed to internal functions. |
A matrix where rows represent metacells and columns represent genes.
data(example_matrix) meta_cells_matrix <- meta_cells( example_matrix ) dim(meta_cells_matrix) meta_cells_matrix[1:6, 1:6]data(example_matrix) meta_cells_matrix <- meta_cells( example_matrix ) dim(meta_cells_matrix) meta_cells_matrix[1:6, 1:6]
Format network table
network_format( network_table, regulators = NULL, targets = NULL, abs_weight = TRUE )network_format( network_table, regulators = NULL, targets = NULL, abs_weight = TRUE )
network_table |
The weight data table of network. |
regulators |
Regulators list. |
targets |
Targets list. |
abs_weight |
Logical value, default is |
Formated network table
data(example_matrix) network_table <- inferCSN(example_matrix) network_format( network_table, regulators = "g1" ) network_format( network_table, regulators = "g1", abs_weight = FALSE ) network_format( network_table, targets = "g3" ) network_format( network_table, regulators = c("g1", "g3"), targets = c("g3", "g5") )data(example_matrix) network_table <- inferCSN(example_matrix) network_format( network_table, regulators = "g1" ) network_format( network_table, regulators = "g1", abs_weight = FALSE ) network_format( network_table, targets = "g3" ) network_format( network_table, regulators = c("g1", "g3"), targets = c("g3", "g5") )
Plot coefficients
plot_coefficient( data, style = "continuous", positive_color = "#3d67a2", negative_color = "#c82926", neutral_color = "#cccccc", bar_width = 0.7, text_size = 3, show_values = TRUE )plot_coefficient( data, style = "continuous", positive_color = "#3d67a2", negative_color = "#c82926", neutral_color = "#cccccc", bar_width = 0.7, text_size = 3, show_values = TRUE )
data |
Input data. |
style |
Plotting style: |
positive_color |
Color for positive weights.
Default is |
negative_color |
Color for negative weights.
Default is |
neutral_color |
Color for weights near zero (used in |
bar_width |
Width of the bars.
Default is |
text_size |
Size of the text for weight values.
Default is |
show_values |
Whether to show weight values on bars.
Default is |
A ggplot object
data(example_matrix) network_table <- inferCSN(example_matrix, targets = "g1") plot_coefficient(network_table) plot_coefficient(network_table, style = "binary")data(example_matrix) network_table <- inferCSN(example_matrix, targets = "g1") plot_coefficient(network_table) plot_coefficient(network_table, style = "binary")
Plot coefficients for multiple targets
plot_coefficients(data, targets = NULL, nrow = NULL, ...)plot_coefficients(data, targets = NULL, nrow = NULL, ...)
data |
Input data. |
targets |
Targets to plot. Default is 'NULL'. |
nrow |
Number of rows for the plot. Default is 'NULL'. |
... |
Other arguments passed to [plot_coefficient]. |
A list of ggplot objects
data(example_matrix) network_table <- inferCSN( example_matrix, targets = c("g1", "g2", "g3") ) plot_coefficients(network_table, show_values = FALSE) plot_coefficients(network_table, targets = "g1")data(example_matrix) network_table <- inferCSN( example_matrix, targets = c("g1", "g2", "g3") ) plot_coefficients(network_table, show_values = FALSE) plot_coefficients(network_table, targets = "g1")
Plot contrast networks
plot_contrast_networks( network_table, degree_value = 0, weight_value = 0, legend_position = "bottom" )plot_contrast_networks( network_table, degree_value = 0, weight_value = 0, legend_position = "bottom" )
network_table |
The weight data table of network. |
degree_value |
Degree value to filter nodes.
Default is |
weight_value |
Weight value to filter edges.
Default is |
legend_position |
The position of legend.
Default is |
A ggplot2 object.
data(example_matrix) network_table <- inferCSN(example_matrix) plot_contrast_networks(network_table[1:50, ])data(example_matrix) network_table <- inferCSN(example_matrix) plot_contrast_networks(network_table[1:50, ])
Plot dynamic networks
plot_dynamic_networks( network_table, celltypes_order, ntop = 10, title = NULL, theme_type = "theme_void", plot_type = "ggplot", layout = "fruchtermanreingold", nrow = 2, figure_save = FALSE, figure_name = NULL, figure_width = 6, figure_height = 6, seed = 1 )plot_dynamic_networks( network_table, celltypes_order, ntop = 10, title = NULL, theme_type = "theme_void", plot_type = "ggplot", layout = "fruchtermanreingold", nrow = 2, figure_save = FALSE, figure_name = NULL, figure_width = 6, figure_height = 6, seed = 1 )
network_table |
The weight data table of network. |
celltypes_order |
The order of cell types. |
ntop |
The number of top genes to plot.
Default is |
title |
The title of figure.
Default is |
theme_type |
The theme of figure.
Could be |
plot_type |
The type of figure.
Could be |
layout |
The layout of figure.
Could be |
nrow |
The number of rows of figure.
Default is |
figure_save |
Whether to save the figure file.
Default is |
figure_name |
The name of figure file.
Default is |
figure_width |
The width of figure.
Default is |
figure_height |
The height of figure.
Default is |
seed |
The seed random use to plot network.
Default is |
A dynamic figure object.
data(example_matrix) network <- inferCSN(example_matrix)[1:100, ] network$celltype <- c( rep("cluster1", 20), rep("cluster2", 20), rep("cluster3", 20), rep("cluster5", 20), rep("cluster6", 20) ) celltypes_order <- c( "cluster5", "cluster3", "cluster2", "cluster1", "cluster6" ) plot_dynamic_networks( network, celltypes_order = celltypes_order ) plot_dynamic_networks( network, celltypes_order = celltypes_order[1:3] ) plot_dynamic_networks( network, celltypes_order = celltypes_order, plot_type = "ggplotly" )data(example_matrix) network <- inferCSN(example_matrix)[1:100, ] network$celltype <- c( rep("cluster1", 20), rep("cluster2", 20), rep("cluster3", 20), rep("cluster5", 20), rep("cluster6", 20) ) celltypes_order <- c( "cluster5", "cluster3", "cluster2", "cluster1", "cluster6" ) plot_dynamic_networks( network, celltypes_order = celltypes_order ) plot_dynamic_networks( network, celltypes_order = celltypes_order[1:3] ) plot_dynamic_networks( network, celltypes_order = celltypes_order, plot_type = "ggplotly" )
Generates visualizations comparing edges of two networks.
plot_edges_comparison( network_table, ground_truth, color_pattern = list(predicted = "gray", ground_truth = "#bb141a", overlap = "#1966ad", total = "#6C757D") )plot_edges_comparison( network_table, ground_truth, color_pattern = list(predicted = "gray", ground_truth = "#bb141a", overlap = "#1966ad", total = "#6C757D") )
network_table |
A data frame of predicted network structure. |
ground_truth |
A data frame of ground truth network. |
color_pattern |
A list of colors for different categories. |
A patchwork plot object containing network edge comparison and distribution plots
data(example_matrix) data("example_ground_truth") network_table <- inferCSN(example_matrix) plot_edges_comparison( network_table, example_ground_truth )data(example_matrix) data("example_ground_truth") network_table <- inferCSN(example_matrix) plot_edges_comparison( network_table, example_ground_truth )
Plot embedding
plot_embedding( matrix, labels = NULL, method = "pca", colors = RColorBrewer::brewer.pal(length(unique(labels)), "Set1"), seed = 1, point_size = 1, cores = 1 )plot_embedding( matrix, labels = NULL, method = "pca", colors = RColorBrewer::brewer.pal(length(unique(labels)), "Set1"), seed = 1, point_size = 1, cores = 1 )
matrix |
Input matrix. |
labels |
Input labels. |
method |
Method to use for dimensionality reduction. |
colors |
Colors to use for the plot. |
seed |
Seed for the random number generator. |
point_size |
Size of the points. |
cores |
Set the number of threads when setting for uwot::umap and Rtsne::Rtsne. |
An embedding plot
data(example_matrix) samples_use <- 1:200 plot_data <- example_matrix[samples_use, ] labels <- sample( c("A", "B", "C", "D", "E"), nrow(plot_data), replace = TRUE ) plot_embedding( plot_data, labels, method = "pca", point_size = 2 ) plot_embedding( plot_data, labels, method = "tsne", point_size = 2 )data(example_matrix) samples_use <- 1:200 plot_data <- example_matrix[samples_use, ] labels <- sample( c("A", "B", "C", "D", "E"), nrow(plot_data), replace = TRUE ) plot_embedding( plot_data, labels, method = "pca", point_size = 2 ) plot_embedding( plot_data, labels, method = "tsne", point_size = 2 )
Plot histogram
plot_histogram( data, binwidth = 0.01, show_border = FALSE, border_color = "black", alpha = 1, theme = "viridis", theme_begin = 0, theme_end = 0.5, theme_direction = -1, legend_position = "right" )plot_histogram( data, binwidth = 0.01, show_border = FALSE, border_color = "black", alpha = 1, theme = "viridis", theme_begin = 0, theme_end = 0.5, theme_direction = -1, legend_position = "right" )
data |
A numeric vector. |
binwidth |
Width of the bins. |
show_border |
Logical value, whether to show border of the bins. |
border_color |
Color of the border. |
alpha |
Alpha value of the bins. |
theme |
Theme of the bins. |
theme_begin |
Begin value of the theme. |
theme_end |
End value of the theme. |
theme_direction |
Direction of the theme. |
legend_position |
The position of legend. |
A ggplot object
data(example_matrix) network_table <- inferCSN(example_matrix) plot_histogram(network_table[, 3])data(example_matrix) network_table <- inferCSN(example_matrix) plot_histogram(network_table[, 3])
Plot network heatmap
plot_network_heatmap( network_table, regulators = NULL, targets = NULL, switch_matrix = TRUE, show_names = FALSE, heatmap_size_lock = TRUE, heatmap_size = 5, heatmap_height = NULL, heatmap_width = NULL, heatmap_title = NULL, heatmap_color = c("#1966ad", "white", "#bb141a"), border_color = "gray", rect_color = NA, anno_width = 1, anno_height = 1, row_anno_type = c("boxplot", "barplot", "histogram", "density", "lines", "points", "horizon"), column_anno_type = c("boxplot", "barplot", "histogram", "density", "lines", "points"), legend_name = "Weight", row_title = "Regulators" )plot_network_heatmap( network_table, regulators = NULL, targets = NULL, switch_matrix = TRUE, show_names = FALSE, heatmap_size_lock = TRUE, heatmap_size = 5, heatmap_height = NULL, heatmap_width = NULL, heatmap_title = NULL, heatmap_color = c("#1966ad", "white", "#bb141a"), border_color = "gray", rect_color = NA, anno_width = 1, anno_height = 1, row_anno_type = c("boxplot", "barplot", "histogram", "density", "lines", "points", "horizon"), column_anno_type = c("boxplot", "barplot", "histogram", "density", "lines", "points"), legend_name = "Weight", row_title = "Regulators" )
network_table |
The weight data table of network. |
regulators |
Regulators list. |
targets |
Targets list. |
switch_matrix |
Whether to weight data table to matrix.
Default is |
show_names |
Whether to show names of row and column.
Default is |
heatmap_size_lock |
Lock the size of heatmap. |
heatmap_size |
The size of heatmap.
Default is |
heatmap_height |
The height of heatmap. |
heatmap_width |
The width of heatmap. |
heatmap_title |
The title of heatmap. |
heatmap_color |
Colors of heatmap. |
border_color |
Default is |
rect_color |
Default is |
anno_width |
Width of annotation. |
anno_height |
Height of annotation. |
row_anno_type |
Default is |
column_anno_type |
Default is |
legend_name |
The name of legend. |
row_title |
The title of row. |
A heatmap
data(example_matrix) data("example_ground_truth") network_table <- inferCSN(example_matrix) p1 <- plot_network_heatmap( example_ground_truth[, 1:3], heatmap_title = "Ground truth", legend_name = "Ground truth" ) p2 <- plot_network_heatmap( network_table, heatmap_title = "inferCSN", legend_name = "inferCSN" ) ComplexHeatmap::draw(p1 + p2) plot_network_heatmap( network_table, show_names = TRUE, rect_color = "gray90", row_anno_type = "density", column_anno_type = "barplot" ) plot_network_heatmap( network_table, regulators = c("g1", "g3", "g5"), targets = c("g3", "g6", "g9"), show_names = TRUE )data(example_matrix) data("example_ground_truth") network_table <- inferCSN(example_matrix) p1 <- plot_network_heatmap( example_ground_truth[, 1:3], heatmap_title = "Ground truth", legend_name = "Ground truth" ) p2 <- plot_network_heatmap( network_table, heatmap_title = "inferCSN", legend_name = "inferCSN" ) ComplexHeatmap::draw(p1 + p2) plot_network_heatmap( network_table, show_names = TRUE, rect_color = "gray90", row_anno_type = "density", column_anno_type = "barplot" ) plot_network_heatmap( network_table, regulators = c("g1", "g3", "g5"), targets = c("g3", "g6", "g9"), show_names = TRUE )
Plot expression data in a scatter plot
plot_scatter( data, smoothing_method = "lm", group_colors = RColorBrewer::brewer.pal(9, "Set1"), title_color = "black", title = NULL, col_title = NULL, row_title = NULL, legend_title = NULL, legend_position = "bottom", margins = "both", marginal_type = NULL, margins_size = 10, compute_correlation = TRUE, compute_correlation_method = "pearson", keep_aspect_ratio = TRUE, facet = FALSE, se = FALSE, pointdensity = TRUE )plot_scatter( data, smoothing_method = "lm", group_colors = RColorBrewer::brewer.pal(9, "Set1"), title_color = "black", title = NULL, col_title = NULL, row_title = NULL, legend_title = NULL, legend_position = "bottom", margins = "both", marginal_type = NULL, margins_size = 10, compute_correlation = TRUE, compute_correlation_method = "pearson", keep_aspect_ratio = TRUE, facet = FALSE, se = FALSE, pointdensity = TRUE )
data |
Input data. |
smoothing_method |
Method for smoothing curve, |
group_colors |
Colors for different groups. |
title_color |
Color for the title. |
title |
Main title for the plot. |
col_title |
Title for the x-axis. |
row_title |
Title for the y-axis. |
legend_title |
Title for the legend. |
legend_position |
The position of legend. |
margins |
The position of marginal figure ("both", "x", "y"). |
marginal_type |
The type of marginal figure ( |
margins_size |
The size of marginal figure, note the bigger size the smaller figure. |
compute_correlation |
Whether to compute and print correlation on the figure. |
compute_correlation_method |
Method to compute correlation ( |
keep_aspect_ratio |
Logical value, whether to set aspect ratio to 1:1. |
facet |
Faceting variable. If setting TRUE, all settings about margins will be inalidation. |
se |
Display confidence interval around smooth. |
pointdensity |
Plot point density when only provide 1 cluster. |
ggplot object
data(example_matrix) test_data <- data.frame( example_matrix[1:200, c(1, 7)], c = c( rep("c1", 40), rep("c2", 40), rep("c3", 40), rep("c4", 40), rep("c5", 40) ) ) p1 <- plot_scatter( test_data ) p2 <- plot_scatter( test_data, marginal_type = "boxplot" ) p1 + p2 p3 <- plot_scatter( test_data, facet = TRUE ) p3 p4 <- plot_scatter( test_data[, 1:2], marginal_type = "histogram" ) p4data(example_matrix) test_data <- data.frame( example_matrix[1:200, c(1, 7)], c = c( rep("c1", 40), rep("c2", 40), rep("c3", 40), rep("c4", 40), rep("c5", 40) ) ) p1 <- plot_scatter( test_data ) p2 <- plot_scatter( test_data, marginal_type = "boxplot" ) p1 + p2 p3 <- plot_scatter( test_data, facet = TRUE ) p3 p4 <- plot_scatter( test_data[, 1:2], marginal_type = "histogram" ) p4
Plot dynamic networks
plot_static_networks( network_table, regulators = NULL, targets = NULL, legend_position = "right" )plot_static_networks( network_table, regulators = NULL, targets = NULL, legend_position = "right" )
network_table |
The weight data table of network. |
regulators |
Regulators list. |
targets |
Targets list. |
legend_position |
The position of legend. |
A ggplot2 object
data(example_matrix) network_table <- inferCSN(example_matrix) plot_static_networks( network_table, regulators = "g1" ) plot_static_networks( network_table, targets = "g1" ) plot_static_networks( network_table, regulators = "g1", targets = "g2" )data(example_matrix) network_table <- inferCSN(example_matrix) plot_static_networks( network_table, regulators = "g1" ) plot_static_networks( network_table, targets = "g1" ) plot_static_networks( network_table, regulators = "g1", targets = "g2" )
Print logo
## S3 method for class 'infercsn_logo' print(x, ...)## S3 method for class 'infercsn_logo' print(x, ...)
x |
Input information. |
... |
Other parameters. |
Print the ASCII logo
Construct network for single target gene
single_network( matrix, regulators, target, cross_validation = FALSE, seed = 1, penalty = "L0", r_squared_threshold = 0, n_folds = 5, verbose = TRUE, ... )single_network( matrix, regulators, target, cross_validation = FALSE, seed = 1, penalty = "L0", r_squared_threshold = 0, n_folds = 5, verbose = TRUE, ... )
matrix |
An expression matrix. |
regulators |
The regulator genes for which to infer the regulatory network. |
target |
The target gene. |
cross_validation |
Whether to use cross-validation.
Default is |
seed |
The random seed for cross-validation.
Default is |
penalty |
The type of regularization.
This can take either one of the following choices: |
r_squared_threshold |
Threshold of R² coefficient.
Default is |
n_folds |
The number of folds for cross-validation.
Default is |
verbose |
Whether to print progress messages.
Default is |
... |
Parameters for other methods. |
A data frame of the single target gene network. The data frame has three columns: regulator, target, and weight.
data(example_matrix) head( single_network( example_matrix, regulators = colnames(example_matrix), target = "g1" ) ) head( single_network( example_matrix, regulators = colnames(example_matrix), target = "g1", cross_validation = TRUE ) ) single_network( example_matrix, regulators = c("g1", "g2", "g3"), target = "g1" ) single_network( example_matrix, regulators = c("g1", "g2"), target = "g1" )data(example_matrix) head( single_network( example_matrix, regulators = colnames(example_matrix), target = "g1" ) ) head( single_network( example_matrix, regulators = colnames(example_matrix), target = "g1", cross_validation = TRUE ) ) single_network( example_matrix, regulators = c("g1", "g2", "g3"), target = "g1" ) single_network( example_matrix, regulators = c("g1", "g2"), target = "g1" )
This function subsamples a matrix using either random sampling or meta cells method.
subsampling( matrix, subsampling_method = c("sample", "meta_cells", "pseudobulk"), subsampling_ratio = 1, seed = 1, verbose = TRUE, ... )subsampling( matrix, subsampling_method = c("sample", "meta_cells", "pseudobulk"), subsampling_ratio = 1, seed = 1, verbose = TRUE, ... )
matrix |
The input matrix to be subsampled. |
subsampling_method |
The method to use for subsampling.
Options are |
subsampling_ratio |
The percent of all samples used for fit_srm.
Default is |
seed |
The random seed for cross-validation.
Default is |
verbose |
Whether to print progress messages.
Default is |
... |
Parameters for other methods. |
The subsampled matrix.
data(example_matrix) data("example_ground_truth") subsample_matrix <- subsampling( example_matrix, subsampling_ratio = 0.5 ) subsample_matrix_2 <- subsampling( example_matrix, subsampling_method = "meta_cells", subsampling_ratio = 0.5, fast_pca = FALSE ) subsample_matrix_3 <- subsampling( example_matrix, subsampling_method = "pseudobulk", subsampling_ratio = 0.5 ) calculate_metrics( inferCSN(example_matrix), example_ground_truth, return_plot = TRUE ) calculate_metrics( inferCSN(subsample_matrix), example_ground_truth, return_plot = TRUE ) calculate_metrics( inferCSN(subsample_matrix_2), example_ground_truth, return_plot = TRUE ) calculate_metrics( inferCSN(subsample_matrix_3), example_ground_truth, return_plot = TRUE )data(example_matrix) data("example_ground_truth") subsample_matrix <- subsampling( example_matrix, subsampling_ratio = 0.5 ) subsample_matrix_2 <- subsampling( example_matrix, subsampling_method = "meta_cells", subsampling_ratio = 0.5, fast_pca = FALSE ) subsample_matrix_3 <- subsampling( example_matrix, subsampling_method = "pseudobulk", subsampling_ratio = 0.5 ) calculate_metrics( inferCSN(example_matrix), example_ground_truth, return_plot = TRUE ) calculate_metrics( inferCSN(subsample_matrix), example_ground_truth, return_plot = TRUE ) calculate_metrics( inferCSN(subsample_matrix_2), example_ground_truth, return_plot = TRUE ) calculate_metrics( inferCSN(subsample_matrix_3), example_ground_truth, return_plot = TRUE )