inferring Cell-Specific gene regulatory Network
inferCSN(
object,
penalty = "L0",
algorithm = "CD",
cross_validation = FALSE,
seed = 1,
n_folds = 10,
subsampling_method = "sample",
subsampling_ratio = 1,
r_threshold = 0,
regulators = NULL,
targets = NULL,
cores = 1,
verbose = TRUE,
...
)
# S4 method for class 'matrix'
inferCSN(
object,
penalty = "L0",
algorithm = "CD",
cross_validation = FALSE,
seed = 1,
n_folds = 10,
subsampling_method = "sample",
subsampling_ratio = 1,
r_threshold = 0,
regulators = NULL,
targets = NULL,
cores = 1,
verbose = TRUE,
...
)
# S4 method for class 'sparseMatrix'
inferCSN(
object,
penalty = "L0",
algorithm = "CD",
cross_validation = FALSE,
seed = 1,
n_folds = 10,
subsampling_method = "sample",
subsampling_ratio = 1,
r_threshold = 0,
regulators = NULL,
targets = NULL,
cores = 1,
verbose = TRUE,
...
)
# S4 method for class 'data.frame'
inferCSN(
object,
penalty = "L0",
algorithm = "CD",
cross_validation = FALSE,
seed = 1,
n_folds = 10,
subsampling_method = "sample",
subsampling_ratio = 1,
r_threshold = 0,
regulators = NULL,
targets = NULL,
cores = 1,
verbose = TRUE,
...
)
The input data for inferCSN
.
The type of regularization, default is L0
.
This can take either one of the following choices: L0
, L0L1
, and L0L2
.
For high-dimensional and sparse data, L0L2
is more effective.
The type of algorithm used to minimize the objective function, default is CD
.
Currently CD
and CDPSI
are supported.
The CDPSI
algorithm may yield better results, but it also increases running time.
Logical value, default is FALSE
, whether to use cross-validation.
The random seed for cross-validation, default is 1
.
The number of folds for cross-validation, default is 10
.
The method to use for subsampling. Options are "sample" or "meta_cells".
The percent of all samples used for sparse_regression
, default is 1
.
Threshold of \(R^2\) or correlation coefficient, default is 0
.
The regulator genes for which to infer the regulatory network.
The target genes for which to infer the regulatory network. 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.
The number of cores to use for parallelization with foreach
, default is 1
.
Logical value, default is TRUE
, whether to print progress messages.
Parameters for other methods.
A data table of regulator-target regulatory relationships
data("example_matrix")
data("example_ground_truth")
network_table_1 <- inferCSN(
example_matrix
)
#> ✔ Running for <dense matrix>.
#> ✔ Checking input parameters.
#> ✔ Using L0 sparse regression model.
#> ✔ Using 1 core.
#> ✔ Run done.
network_table_2 <- inferCSN(
example_matrix,
cores = 2
)
#> ✔ Running for <dense matrix>.
#> ✔ Checking input parameters.
#> ✔ Using L0 sparse regression model.
#> ✔ Using 2 cores.
#> ✔ Run done.
head(network_table_1)
#> regulator target weight
#> 1 g17 g18 0.5443805
#> 2 g18 g1 -0.5390453
#> 3 g16 g15 0.4818454
#> 4 g17 g16 0.4599054
#> 5 g15 g14 0.4516497
#> 6 g14 g13 0.4483912
identical(
network_table_1,
network_table_2
)
#> [1] TRUE
inferCSN(
example_matrix,
regulators = c("g1", "g2"),
targets = c("g3", "g4")
)
#> ✔ Running for <dense matrix>.
#> ✔ Checking input parameters.
#> ✔ Using 2 regulator(s).
#> ✔ Using 2 target(s).
#> ✔ Using L0 sparse regression model.
#> ✔ Using 1 core.
#> ✔ Run done.
#> regulator target weight
#> 1 g2 g3 0.8504059
#> 2 g2 g4 0.7058242
#> 3 g1 g4 -0.2941758
#> 4 g1 g3 -0.1495941
inferCSN(
example_matrix,
regulators = c("g1", "g2"),
targets = c("g3", "g0")
)
#> ✔ Running for <dense matrix>.
#> ✔ Checking input parameters.
#> ✔ Using 2 regulator(s).
#> ! Warning: 1 out of 2 candidate targets are in the input matrix.
#> ✔ Using L0 sparse regression model.
#> ✔ Using 1 core.
#> ✔ Run done.
#> regulator target weight
#> 1 g2 g3 0.8504059
#> 2 g1 g3 -0.1495941
if (FALSE) { # \dontrun{
network_table_07 <- inferCSN(
example_matrix,
r_threshold = 0.7
)
calculate_metrics(
network_table_1,
example_ground_truth,
plot = TRUE
)
calculate_metrics(
network_table_07,
example_ground_truth,
plot = TRUE
)
} # }
if (FALSE) { # \dontrun{
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
)
} # }