Title: | Single-Cell RNA Sequencing Pseudo-Time Based Gene Regulatory Network Inference |
---|---|
Description: | Inference and visualize gene regulatory network based on single-cell RNA sequencing pseudo-time information. |
Authors: | Gangcai Xie |
Maintainer: | Gangcai Xie <[email protected]> |
License: | MIT + file LICENSE |
Version: | 0.3.5 |
Built: | 2025-02-18 05:43:19 UTC |
Source: | https://github.com/cran/pGRN |
Based on single-cell pseudotime information, get the sliding window average expression, and then standard normlize the expression for each gene
data_transform(data, pseudotime, slide_window_size = 100, slide_step_size = 50)
data_transform(data, pseudotime, slide_window_size = 100, slide_step_size = 50)
data |
expression matrix data |
pseudotime |
list of pseudotime |
slide_window_size |
sliding window size |
slide_step_size |
sliding window step size |
Transformed new matrix
data <- matrix(1,100,1000) ptime <- seq(1:1000) data_transform(data, ptime, slide_window_size=100, slide_step_size=50)
data <- matrix(1,100,1000) ptime <- seq(1:1000) data_transform(data, ptime, slide_window_size=100, slide_step_size=50)
Get bidirectional DTW distance.
get_dtw_dist_bidirectional(x, y)
get_dtw_dist_bidirectional(x, y)
x |
list of x input |
y |
list of y input |
numeric
get_dtw_dist_bidirectional(c(1:1000),c(1:1000))
get_dtw_dist_bidirectional(c(1:1000),c(1:1000))
Get DTW distance matrix for all genes using pseudotime based sliding window transfromation, parallel computing allowed.
get_dtw_dist_mat( data, ptime, slide_window_size = 50, slide_step_size = 25, cores = 2 )
get_dtw_dist_mat( data, ptime, slide_window_size = 50, slide_step_size = 25, cores = 2 )
data |
gene expression matrix (Gene * Cells) |
ptime |
pseudotime matched with the column cells of the gene expression matrix |
slide_window_size |
sliding window size |
slide_step_size |
sliding window step size |
cores |
number of cores for parallel computing |
bidirectional DTW distance matrix
example_data <- pGRNDB expression_matrix <- example_data[["expression"]] pseudotime_list <- example_data[["ptime"]]$PseudoTime dtw_dist_matrix <- get_dtw_dist_mat(expression_matrix, pseudotime_list, cores=1)
example_data <- pGRNDB expression_matrix <- example_data[["expression"]] pseudotime_list <- example_data[["ptime"]]$PseudoTime dtw_dist_matrix <- get_dtw_dist_mat(expression_matrix, pseudotime_list, cores=1)
Get sub-networks based on given adjacency data.frame input
get_networks( data, centrality_degree_mod = "out", components_mod = "weak", network_min_genes = 10 )
get_networks( data, centrality_degree_mod = "out", components_mod = "weak", network_min_genes = 10 )
data |
adjacency data.frame |
centrality_degree_mod |
mode of centrality degree for popularity calculation |
components_mod |
mode of sub-network extraction methods |
network_min_genes |
minimal number of gene elements required for extracted sub-networks |
list of tabl_graph objects
example_data <- pGRNDB expression_matrix <- example_data[["expression"]] pseudotime_list <- example_data[["ptime"]]$PseudoTime dtw_dist_matrix <- get_dtw_dist_mat(expression_matrix, pseudotime_list, cores=1) adj_df <- matrix2adj(dtw_dist_matrix) get_networks(adj_df,network_min_genes=5)
example_data <- pGRNDB expression_matrix <- example_data[["expression"]] pseudotime_list <- example_data[["ptime"]]$PseudoTime dtw_dist_matrix <- get_dtw_dist_mat(expression_matrix, pseudotime_list, cores=1) adj_df <- matrix2adj(dtw_dist_matrix) get_networks(adj_df,network_min_genes=5)
Convert distance matrix to adjacency dataframe for network construction.
matrix2adj(data, quantile_cutoff = 5)
matrix2adj(data, quantile_cutoff = 5)
data |
distance matrix |
quantile_cutoff |
an integer value (1-99) for quantile cutoff |
adjacency dataframe (with columns "from, to, distance,direction, similarity")
example_data <- pGRNDB expression_matrix <- example_data[["expression"]] pseudotime_list <- example_data[["ptime"]]$PseudoTime dtw_dist_matrix <- get_dtw_dist_mat(expression_matrix, pseudotime_list, cores=1) adj_df <- matrix2adj(dtw_dist_matrix)
example_data <- pGRNDB expression_matrix <- example_data[["expression"]] pseudotime_list <- example_data[["ptime"]]$PseudoTime dtw_dist_matrix <- get_dtw_dist_mat(expression_matrix, pseudotime_list, cores=1) adj_df <- matrix2adj(dtw_dist_matrix)
Given a distance matrix, calculate gene modules based on hierarchical clustering method and then get module level networks
module_networks( data, k = 10, quantile_cutoff = 10, centrality_degree_mod = "out", components_mod = "weak", network_min_genes = 10 )
module_networks( data, k = 10, quantile_cutoff = 10, centrality_degree_mod = "out", components_mod = "weak", network_min_genes = 10 )
data |
distance matrix |
k |
number of gene clusters for module inference |
quantile_cutoff |
distance cutoff based on quantile(1-99) for edge identification |
centrality_degree_mod |
"in" or "out" for nodes popularity calculation |
components_mod |
"weak" or "strong" for sub-network components inference |
network_min_genes |
minial number of genes required for a network |
a list networks for each module
example_data <- pGRNDB expression_matrix <- example_data[["expression"]] pseudotime_list <- example_data[["ptime"]]$PseudoTime dtw_dist_matrix <- get_dtw_dist_mat(expression_matrix, pseudotime_list, cores=1) nets <- module_networks(dtw_dist_matrix,k=1,quantile_cutoff=50) plot_network(nets[["module1"]])
example_data <- pGRNDB expression_matrix <- example_data[["expression"]] pseudotime_list <- example_data[["ptime"]]$PseudoTime dtw_dist_matrix <- get_dtw_dist_mat(expression_matrix, pseudotime_list, cores=1) nets <- module_networks(dtw_dist_matrix,k=1,quantile_cutoff=50) plot_network(nets[["module1"]])
Given single cell matrix and pseudotime, construct gene regulatory network (GRN)
pGRN( expression_matrix, pseudotime_list, method = "DTW", slide_window_size = 20, slide_step_size = 10, centrality_degree_mod = "out", components_mod = "weak", network_min_genes = 10, quantile_cutoff = 5, order = 1, cores = 1 )
pGRN( expression_matrix, pseudotime_list, method = "DTW", slide_window_size = 20, slide_step_size = 10, centrality_degree_mod = "out", components_mod = "weak", network_min_genes = 10, quantile_cutoff = 5, order = 1, cores = 1 )
expression_matrix |
expression matrix data |
pseudotime_list |
list of pseudotime |
method |
method for GRN construction: DTW, granger |
slide_window_size |
sliding window size |
slide_step_size |
sliding window step size |
centrality_degree_mod |
(for DTW method) mode of centrality degree for popularity calculation |
components_mod |
(for DTW method) mode of sub-network extraction methods (weak or strong) |
network_min_genes |
minimal number of gene elements required for extracted sub-networks |
quantile_cutoff |
an integer value (1-99) for quantile cutoff |
order |
(for granger method) integer specifying the order of lags to include in the auxiliary regression |
cores |
number of cores for parallel computing |
a list of tabl_graph objects
example_data <- pGRNDB expression_matrix <- example_data[["expression"]] pseudotime_list <- example_data[["ptime"]]$PseudoTime # try DTW method nets <- pGRN(expression_matrix, pseudotime_list, method= "DTW", quantile_cutoff=50, cores=1) plot_network(nets[[1]]) # plot the network interactively plot_network_i(nets[[1]])
example_data <- pGRNDB expression_matrix <- example_data[["expression"]] pseudotime_list <- example_data[["ptime"]]$PseudoTime # try DTW method nets <- pGRN(expression_matrix, pseudotime_list, method= "DTW", quantile_cutoff=50, cores=1) plot_network(nets[[1]]) # plot the network interactively plot_network_i(nets[[1]])
A list with expression dataframe and pseudotime dataframe
pGRNDB
pGRNDB
pGRNDB
A list with items expression and ptime
data frame of single cell expression
pseudotime of the single cells
...
pGRN
Plot stationary network through ggraph
plot_network(graph, ...)
plot_network(graph, ...)
graph |
a tbl_graph object |
... |
other parameters for ggraph |
ggraph
example_data <- pGRNDB expression_matrix <- example_data[["expression"]] pseudotime_list <- example_data[["ptime"]]$PseudoTime dtw_dist_matrix <- get_dtw_dist_mat(expression_matrix, pseudotime_list, cores=1) nets <- module_networks(dtw_dist_matrix,k=1,quantile_cutoff=50) plot_network(nets[["module1"]])
example_data <- pGRNDB expression_matrix <- example_data[["expression"]] pseudotime_list <- example_data[["ptime"]]$PseudoTime dtw_dist_matrix <- get_dtw_dist_mat(expression_matrix, pseudotime_list, cores=1) nets <- module_networks(dtw_dist_matrix,k=1,quantile_cutoff=50) plot_network(nets[["module1"]])
Plot interactive network based on igraph layout input
plot_network_i(graph, save_file = NULL)
plot_network_i(graph, save_file = NULL)
graph |
igraph layout object |
save_file |
file name of the saved file, not save if NULL |
visNetwork htmlwidget
example_data <- pGRNDB expression_matrix <- example_data[["expression"]] pseudotime_list <- example_data[["ptime"]]$PseudoTime dtw_dist_matrix <- get_dtw_dist_mat(expression_matrix, pseudotime_list, cores=1) nets <- module_networks(dtw_dist_matrix,k=1,quantile_cutoff=50) plot_network_i(nets[["module1"]])
example_data <- pGRNDB expression_matrix <- example_data[["expression"]] pseudotime_list <- example_data[["ptime"]]$PseudoTime dtw_dist_matrix <- get_dtw_dist_mat(expression_matrix, pseudotime_list, cores=1) nets <- module_networks(dtw_dist_matrix,k=1,quantile_cutoff=50) plot_network_i(nets[["module1"]])
Use DTW to calcuate gene-gene distance based on their expression and pseudotime
run_dtw( expression_matrix, pseudotime_list, slide_window_size = 50, slide_step_size = 25, quantile_cutoff = 5, cores = 1 )
run_dtw( expression_matrix, pseudotime_list, slide_window_size = 50, slide_step_size = 25, quantile_cutoff = 5, cores = 1 )
expression_matrix |
expression matrix data |
pseudotime_list |
list of pseudotime |
slide_window_size |
sliding window size |
slide_step_size |
sliding window step size |
quantile_cutoff |
an integer value (1-99) for quantile cutoff |
cores |
number of cores for parallel computing |
adjacency dataframe (with columns "from, to, distance,direction, similarity")
example_data <- pGRNDB expression_matrix <- example_data[["expression"]] pseudotime_list <- example_data[["ptime"]]$PseudoTime adj_df <- run_dtw(expression_matrix, pseudotime_list, quantile_cutoff=50, cores=1)
example_data <- pGRNDB expression_matrix <- example_data[["expression"]] pseudotime_list <- example_data[["ptime"]]$PseudoTime adj_df <- run_dtw(expression_matrix, pseudotime_list, quantile_cutoff=50, cores=1)
Based on single-cell gene expression matrix and pseudotime, calculate Granger-causality Test based gene-gene regulatory relationship
run_granger_test( data, ptime, slide_window_size = 20, slide_step_size = 10, pvalue_cutoff = 0.01, order = 1, ... )
run_granger_test( data, ptime, slide_window_size = 20, slide_step_size = 10, pvalue_cutoff = 0.01, order = 1, ... )
data |
gene expression matrix (Gene * Cells) |
ptime |
pseudotime matched with the column cells of the gene expression matrix |
slide_window_size |
sliding window size |
slide_step_size |
sliding window step size |
pvalue_cutoff |
cutoff for the pvalue from transfer entropy test |
order |
integer specifying the order of lags to include in the auxiliary regression |
... |
other parameters for grangertest function in lmtest |
adjacency data frame
example_data <- pGRNDB expression_matrix <- example_data[["expression"]] pseudotime_list <- example_data[["ptime"]]$PseudoTime gt_adj_df <- run_granger_test(expression_matrix, pseudotime_list)
example_data <- pGRNDB expression_matrix <- example_data[["expression"]] pseudotime_list <- example_data[["ptime"]]$PseudoTime gt_adj_df <- run_granger_test(expression_matrix, pseudotime_list)
Get sliding windows average values for given vector/list
slideWindows(data, window = 2, step = 1)
slideWindows(data, window = 2, step = 1)
data |
list of expression |
window |
sliding window size |
step |
sliding window step size |
list/vector of sliding windows with average expression value
slideWindows(c(1:1000),window=200,step=100) slideWindows(c(1:1000),window=100,step=50)
slideWindows(c(1:1000),window=200,step=100) slideWindows(c(1:1000),window=100,step=50)