inferCSN | R Documentation |
inferring Cell-Specific gene regulatory Network
inferCSN(
object,
penalty = "L0",
algorithm = "CD",
cross_validation = FALSE,
seed = 1,
n_folds = 10,
percent_samples = 1,
r_threshold = 0,
regulators = NULL,
targets = NULL,
regulators_num = NULL,
cores = 1,
verbose = TRUE,
...
)
## S4 method for signature 'matrix'
inferCSN(
object,
penalty = "L0",
algorithm = "CD",
cross_validation = FALSE,
seed = 1,
n_folds = 10,
percent_samples = 1,
r_threshold = 0,
regulators = NULL,
targets = NULL,
regulators_num = NULL,
cores = 1,
verbose = TRUE,
...
)
## S4 method for signature 'sparseMatrix'
inferCSN(
object,
penalty = "L0",
algorithm = "CD",
cross_validation = FALSE,
seed = 1,
n_folds = 10,
percent_samples = 1,
r_threshold = 0,
regulators = NULL,
targets = NULL,
regulators_num = NULL,
cores = 1,
verbose = TRUE,
...
)
## S4 method for signature 'data.frame'
inferCSN(
object,
penalty = "L0",
algorithm = "CD",
cross_validation = FALSE,
seed = 1,
n_folds = 10,
percent_samples = 1,
r_threshold = 0,
regulators = NULL,
targets = NULL,
regulators_num = NULL,
cores = 1,
verbose = TRUE,
...
)
object |
The input data for |
penalty |
The type of regularization, default is |
algorithm |
The type of algorithm used to minimize the objective function, default is |
cross_validation |
Logical value, default is |
seed |
The random seed for cross-validation, default is |
n_folds |
The number of folds for cross-validation, default is |
percent_samples |
The percent of all samples used for |
r_threshold |
Threshold of |
regulators |
The regulator genes for which to infer the regulatory network. |
targets |
The target genes for which to infer the regulatory network. |
regulators_num |
The number of non-zore coefficients, this value will affect the final performance. The maximum support size at which to terminate the regularization path. Recommend setting this to a small fraction of min(n,p) (e.g. 0.05 * min(n,p)) as L0 regularization typically selects a small portion of non-zeros. |
cores |
The number of cores to use for parallelization with |
verbose |
Logical value, default is |
... |
Parameters for other methods. |
A data table of regulator-target regulatory relationships
data("example_matrix")
network_table_1 <- inferCSN(
example_matrix
)
head(network_table_1)
network_table_2 <- inferCSN(
example_matrix,
cores = 2
)
identical(
network_table_1,
network_table_2
)
inferCSN(
example_matrix,
regulators = c("g1", "g2"),
targets = c("g3", "g4")
)
inferCSN(
example_matrix,
regulators = c("g1", "g2"),
targets = c("g3", "g0")
)
inferCSN(
example_matrix,
regulators = c("g1", "g0"),
targets = c("g2", "g3")
)
inferCSN(
example_matrix,
regulators = c("g1"),
targets = c("g2")
)
## Not run:
data("example_matrix")
network_table <- inferCSN(example_matrix)
head(network_table)
network_table_sparse_1 <- inferCSN(
as(example_matrix, "sparseMatrix")
)
head(network_table_sparse_1)
network_table_sparse_2 <- inferCSN(
as(example_matrix, "sparseMatrix"),
cores = 2
)
identical(
network_table,
network_table_sparse_1
)
identical(
network_table_sparse_1,
network_table_sparse_2
)
plot_scatter(
data.frame(
network_table$weight,
network_table_sparse_1$weight
),
legend_position = "none"
)
plot_weight_distribution(
network_table
) + plot_weight_distribution(
network_table_sparse_1
)
## End(Not run)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.