Title: | Molecular Signature identification using Biclustering |
---|---|
Description: | This package is a implementation of biclustering ensemble method MoSBi (Molecular signature Identification from Biclustering). MoSBi provides standardized interfaces for biclustering results and can combine their results with a multi-algorithm ensemble approach to compute robust ensemble biclusters on molecular omics data. This is done by computing similarity networks of biclusters and filtering for overlaps using a custom error model. After that, the louvain modularity it used to extract bicluster communities from the similarity network, which can then be converted to ensemble biclusters. Additionally, MoSBi includes several network visualization methods to give an intuitive and scalable overview of the results. MoSBi comes with several biclustering algorithms, but can be easily extended to new biclustering algorithms. |
Authors: | Tim Daniel Rose [cre, aut], Josch Konstantin Pauling [aut], Nikolai Koehler [aut] |
Maintainer: | Tim Daniel Rose <[email protected]> |
License: | AGPL-3 + file LICENSE |
Version: | 1.13.0 |
Built: | 2024-10-30 08:55:27 UTC |
Source: | https://github.com/bioc/mosbi |
Can be used for .g. histograms.
algohistogram(bic)
algohistogram(bic)
bic |
A list of bicluster objects. |
A character vector with the extracted biclustering algorithms used for each bicluster of the input list.
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # algohistogram(bics)
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # algohistogram(bics)
All values lower than the threshold will be replaced by a 0.
apply_threshold(bic_net)
apply_threshold(bic_net)
bic_net |
An object of class |
An adjacency matrix with the applied threshold.
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # bn <- bicluster_network(bics, m) # apply_threshold(bn)
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # bn <- bicluster_network(bics, m) # apply_threshold(bn)
All values lower than the threshold will be replaced by a 0.
## S4 method for signature 'bicluster_net' apply_threshold(bic_net)
## S4 method for signature 'bicluster_net' apply_threshold(bic_net)
bic_net |
An object of class |
An adjacency matrix with the applied threshold.
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # bn <- bicluster_network(bics, m) # apply_threshold(bn)
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # bn <- bicluster_network(bics, m) # apply_threshold(bn)
All values lower than the threshold will be replaced by a 0.
## S4 method for signature 'cooccurrence_net' apply_threshold(bic_net)
## S4 method for signature 'cooccurrence_net' apply_threshold(bic_net)
bic_net |
An object of class |
An adjacency matrix with the applied threshold.
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # fn <- feature_network(bics, m) # apply_threshold(fn)
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # fn <- feature_network(bics, m) # apply_threshold(fn)
Given a list of bicluster objects (bicluster
),
the function counts the occurance of all elements in the biclusters.
attr_overlap(bics, named)
attr_overlap(bics, named)
bics |
A list of |
named |
Boolean, indicating, if all bicluster objects have names. |
A Data Frame with the counts oof all elements.
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # attr_overlap(bics, named=FALSE)
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # attr_overlap(bics, named=FALSE)
The function generates co-occurance networks for all the attributes.
E.g. if MARGIN="column"
, for each column, a oc-occurance matrix
of rows is generated, which includes all biclusters, where the
column element is present.
attribute_graph(bics, m, MARGIN = "column")
attribute_graph(bics, m, MARGIN = "column")
bics |
A list of |
m |
The matrix used for biclustering. |
MARGIN |
|
A list of numeric matrices.
If MARGIN="column"
("row"
), the list has a
length of ncol(m)
(nrow(m)
)
and each matrix the dimensions of c(nrow(m),
nrow(m))
(c(ncol(m), ncol(m))
)
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # attribute_graph(bics, m)
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # attribute_graph(bics, m)
For a adjacency matrix as computed by full_graph
,
the function computes how many row-column interactions connect
rows (columns) to columns (rows) of a specific class/category.
attributeConnector(mat, otherclasses, useOther = FALSE)
attributeConnector(mat, otherclasses, useOther = FALSE)
mat |
A adjacency matrix with bipartite interactions as
computed by |
otherclasses |
A logical vector indicating two classes of elements in rows (columns). |
useOther |
Logical indicating if the attributes, that
are classified appear first in the matrix ( |
A DataFrame that holds the total degree of every
attribute (row/column) and the fraction of the degree that
connects only to elements of class True
(from
parameter otherclasses
).
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # fn <- feature_network(bics, m) # attributeConnector(apply_threshold(fn), # otherclasses=c(rep(FALSE, 100), rep(TRUE, 100)))
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # fn <- feature_network(bics, m) # attributeConnector(apply_threshold(fn), # otherclasses=c(rep(FALSE, 100), rep(TRUE, 100)))
Uses the stats::heatmap
function.
bicluster_heatmap(bic, m, ...)
bicluster_heatmap(bic, m, ...)
bic |
A bicluster object. |
m |
The matrix, that was used for the biclustering. (Works only if matrix has row-/colnames.) |
... |
Arguments forwarded to |
A plot object
m <- matrix(c(1,2,3,4), nrow=2) rownames(m) <- c("r1", "r2") rownames(m) <- c("c1", "c2") bicluster_heatmap(bicluster(row=c(1,2), column=c(1,2)), m)
m <- matrix(c(1,2,3,4), nrow=2) rownames(m) <- c("r1", "r2") rownames(m) <- c("c1", "c2") bicluster_heatmap(bicluster(row=c(1,2), column=c(1,2)), m)
Uses the stats::heatmap
function.
## S4 method for signature 'bicluster,matrix' bicluster_heatmap(bic, m, ...)
## S4 method for signature 'bicluster,matrix' bicluster_heatmap(bic, m, ...)
bic |
A bicluster object. |
m |
The matrix, that was used for the biclustering. (Works only if matrix has row-/colnames.) |
... |
Arguments forwarded to |
A plot object
m <- matrix(c(1,2,3,4), nrow=2) rownames(m) <- c("r1", "r2") rownames(m) <- c("c1", "c2") bicluster_heatmap(bicluster(row=c(1,2), column=c(1,2)), m)
m <- matrix(c(1,2,3,4), nrow=2) rownames(m) <- c("r1", "r2") rownames(m) <- c("c1", "c2") bicluster_heatmap(bicluster(row=c(1,2), column=c(1,2)), m)
The function converts a bicluster_net
object into an igraph
graph object.
The threshold
is used as a cutoff for the edges of the network.
bicluster_net_to_igraph(bic_net)
bicluster_net_to_igraph(bic_net)
bic_net |
An object of class |
An igraph::graph
object.
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # bn <- bicluster_network(bics, m) # bicluster_net_to_igraph(bn)
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # bn <- bicluster_network(bics, m) # bicluster_net_to_igraph(bn)
The function converts a bicluster_net
object into an igraph
graph object.
The threshold
is used as a cutoff for the edges of the network.
## S4 method for signature 'bicluster_net' bicluster_net_to_igraph(bic_net)
## S4 method for signature 'bicluster_net' bicluster_net_to_igraph(bic_net)
bic_net |
An object of class |
An igraph::graph
object.
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # bn <- bicluster_network(bics, m) # bicluster_net_to_igraph(bn)
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # bn <- bicluster_network(bics, m) # bicluster_net_to_igraph(bn)
Object that is returned e.g. be the
function bicluster_network
.
adjacency_matrix
Adjacency matrix of bicluster similarities.
threshold
Estimated threshold for the bicluster similarity adjacency
matrix. All values lower than that in the matrix should be discarded.
(Note that the indicated threshold is not applied
to the adjacency_matrix
)
algorithms
List of algorithms that contributed to this bicluster network.
bicluster_net(adjacency_matrix=matrix(seq(1:16), nrow=4), threshold=4)
bicluster_net(adjacency_matrix=matrix(seq(1:16), nrow=4), threshold=4)
The function computes a bicluster network based on a selected similarity metric. A similarity cut-off is calculated using randomized biclusters (the bicluster size distribution is kept).
bicluster_network( bics, mat, n_randomizations = 5, MARGIN = "both", metric = 4, n_steps = 100, plot_edge_dist = TRUE, sn_ratio = TRUE, error_threshold = 0.05, return_plot_data = FALSE, prob_scale = FALSE, prl = FALSE )
bicluster_network( bics, mat, n_randomizations = 5, MARGIN = "both", metric = 4, n_steps = 100, plot_edge_dist = TRUE, sn_ratio = TRUE, error_threshold = 0.05, return_plot_data = FALSE, prob_scale = FALSE, prl = FALSE )
bics |
A list of bicluster objects. |
mat |
The matrix used for biclustering. |
n_randomizations |
The number of randomizations for cut-off estimation. (The mean of all randomizations is used). |
MARGIN |
Margin over which the similarity is computed. Can be "row", "column", "mean" (In this case the mean of row and column similarity is used) or "both" (In this case the similarity between all the datapoints of biclusters is used). |
metric |
The similarity metric same as
in |
n_steps |
Number of points where the difference between randomizations and the real data is evaluated. |
plot_edge_dist |
Show the plots for cut-off estimation with the error model. |
sn_ratio |
If |
error_threshold |
If |
return_plot_data |
Please do not use outside of the package. |
prob_scale |
Scale similarity by the probability of an
overlap equal of higher to the observed one. The scaling is
done by multiplying the similarity
with |
prl |
Compute the similarity matrix using multiple
cores (works only for |
An object of class bicluster_net
.
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # bicluster_network(bics, m)
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # bicluster_network(bics, m)
Convert a bicluster object to an acutal submatrix of the original matrix.
bicluster_to_matrix(m, bic)
bicluster_to_matrix(m, bic)
m |
Matrix on which the bicluster was computed |
bic |
Bicluster object |
A matrix.
bicluster_to_matrix(matrix(seq(1:16), nrow=4), bicluster(row=c(1,2), column=c(1,2)))
bicluster_to_matrix(matrix(seq(1:16), nrow=4), bicluster(row=c(1,2), column=c(1,2)))
Convert a bicluster object to an acutal submatrix of the original matrix.
## S4 method for signature 'matrix,bicluster' bicluster_to_matrix(m, bic)
## S4 method for signature 'matrix,bicluster' bicluster_to_matrix(m, bic)
m |
Matrix on which the bicluster was computed |
bic |
Bicluster object |
A matrix.
#' @examples bicluster_to_matrix(matrix(seq(1:16), nrow=4), bicluster(row=c(1,2), column=c(1,2)))
A S4 class to store biclusters.
row
A vector of row.
column
A vector of columns.
rowname
A vector of names for the rows in row
.
colname
A vector of names for the columns in column
.
algorithm
Algorithm that predicted this bicluster.
bicluster(row=c(1,2), column=c(1,2), rowname=c("a", "b"), colname=c("e", "f"))
bicluster(row=c(1,2), column=c(1,2), rowname=c("a", "b"), colname=c("e", "f"))
Throw an error, if a matrix has not both row- and colnames.
check_names(m)
check_names(m)
m |
A matrix. |
Throws error, if matrix has no row- and column names.
m <- matrix(c(1,2,3,4), nrow=2) rownames(m) <- c("r1", "r2") colnames(m) <- c("c1", "c2") check_names(m)
m <- matrix(c(1,2,3,4), nrow=2) rownames(m) <- c("r1", "r2") colnames(m) <- c("c1", "c2") check_names(m)
Clean a list of biclusters, by returning only the valid ones,
clean_bicluster_list(bics)
clean_bicluster_list(bics)
bics |
A list of bicluster objects. |
A lis tof bicluster objects
b <- list(bicluster(row=c(1,2,3,4), column=c(1,2,3,4)), bicluster(row=c(3,4,5,6), column=c(3,4,5,6))) clean_bicluster_list(b)
b <- list(bicluster(row=c(1,2,3,4), column=c(1,2,3,4)), bicluster(row=c(3,4,5,6), column=c(3,4,5,6))) clean_bicluster_list(b)
Can be used for e.g. histograms.
colhistogram(bic)
colhistogram(bic)
bic |
A list of bicluster objects. |
A vector with the lenghts of the columns in every bicluster object.
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # colhistogram(bics)
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # colhistogram(bics)
The function converts a cooccurrence_net
object into an
igraph graph object.
The threshold
is used as a cutoff for the edges of the network.
cooccurrence_net_to_igraph(occ_net)
cooccurrence_net_to_igraph(occ_net)
occ_net |
An object of class |
An igraph::graph
object.
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # fn <- feature_network(bics, m) # cooccurrence_net_to_igraph(fn)
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # fn <- feature_network(bics, m) # cooccurrence_net_to_igraph(fn)
The function converts a cooccurrence_net
object into an
igraph graph object.
The threshold
is used as a cutoff for the edges of the network.
## S4 method for signature 'cooccurrence_net' cooccurrence_net_to_igraph(occ_net)
## S4 method for signature 'cooccurrence_net' cooccurrence_net_to_igraph(occ_net)
occ_net |
An object of class |
An igraph::graph
object.
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # fn <- feature_network(bics, m) # cooccurrence_net_to_igraph(fn)
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # fn <- feature_network(bics, m) # cooccurrence_net_to_igraph(fn)
Object that is returned e.g. be the function feature_network
.
adjacency_matrix
Adjacency matrix of row- and column-element co-occurences.
threshold
Estimated threshold for the co-occurence adjacency matrix.
All values lower than that in the matrix should be
discarded. (Note that the indicated threshold is not
applied to the adjacency_matrix
)
cooccurrence_net(adjacency_matrix=matrix(seq(1:16), nrow=4), threshold=4)
cooccurrence_net(adjacency_matrix=matrix(seq(1:16), nrow=4), threshold=4)
Subsetting of R matrices within c++.
cpp_matrix_subsetting(m, bic)
cpp_matrix_subsetting(m, bic)
m |
A numeric matrix |
bic |
A bicluster object. |
Matrix subset.
cpp_matrix_subsetting(matrix(seq(1:16), nrow=4), bicluster(row=c(1,2), column=c(1,2)))
cpp_matrix_subsetting(matrix(seq(1:16), nrow=4), bicluster(row=c(1,2), column=c(1,2)))
Finds the highest element in a list of bicluster objects.
detect_elements(bics, MARGIN = "row")
detect_elements(bics, MARGIN = "row")
bics |
A list of bicluster objects. |
MARGIN |
Choose if the distance is computed over |
Return highest row or column index from a list of biclusters.
b <- list(bicluster(row=c(1,2,3,4), column=c(1,2,3,4)), bicluster(row=c(3,4,5,6), column=c(3,4,5,6))) detect_elements(b)
b <- list(bicluster(row=c(1,2,3,4), column=c(1,2,3,4)), bicluster(row=c(3,4,5,6), column=c(3,4,5,6))) detect_elements(b)
Get the dimensions of a bicluster.
## S4 method for signature 'bicluster' dim(x)
## S4 method for signature 'bicluster' dim(x)
x |
A bicluster object. |
A numeric vector with the lengths of the rows and columns of the bicluster.
dim(bicluster(row=c(1,2), column=c(1,2)))
dim(bicluster(row=c(1,2), column=c(1,2)))
This function computes a distance matrix between biclusters using different dissimilarity metrics.
distance_matrix(bics, MARGIN = "row", metric = 1L)
distance_matrix(bics, MARGIN = "row", metric = 1L)
bics |
A list of bicluster objects. |
MARGIN |
Choose if the distance is computed over |
metric |
Integer indicating which metric is used. 1: Bray-Curtis dissimilarity (default), 2: Jaccard distance, 3: 1-overlap coefficient 4: 1 - Fowlkes–Mallows index. |
A numeric matrix of the dissimilarities between all given biclusters.
b <- list(bicluster(row=c(1,2,3,4), column=c(1,2,3,4)), bicluster(row=c(3,4,5,6), column=c(3,4,5,6))) distance_matrix(b)
b <- list(bicluster(row=c(1,2,3,4), column=c(1,2,3,4)), bicluster(row=c(3,4,5,6), column=c(3,4,5,6))) distance_matrix(b)
After calculation of communities with
the get_louvain_communities
function, the result can be
converted into a list of bicluster
objects with this function.
Only biclusters are returned which have a minimum dimension of 2x2.
ensemble_biclusters( coms, bics, mat, row_threshold = 0.1, col_threshold = 0.1, threshold_sorted = FALSE )
ensemble_biclusters( coms, bics, mat, row_threshold = 0.1, col_threshold = 0.1, threshold_sorted = FALSE )
coms |
A list of communities ( |
bics |
The list biclusters that was used for calculation
with |
mat |
The numeric matrix, that was used for biclustering. |
row_threshold |
Minimum fraction of biclusters of a community in which a row needs to occur so that it will be part of the outputted ensemble bicluster. |
col_threshold |
Minimum fraction of biclusters of a community in which a column needs to occur so that it will be part of the outputted ensemble bicluster. |
threshold_sorted |
Return the rows and columns in sorted by decreasing fraction. |
A list of bicluster
objects.
b <- list(bicluster(row=c(1,2,3,4), column=c(1,2,3,4)), bicluster(row=c(3,4,5,6), column=c(3,4,5,6))) # m <- matrix(runif(100), nrow=10) # tm = matrix(c(0,1,1,0), nrow=2) # bn <- list(bicluster_net(adjacency_matrix=tm, threshold=.5)) # ensemble_biclusters(bn, b, m)
b <- list(bicluster(row=c(1,2,3,4), column=c(1,2,3,4)), bicluster(row=c(3,4,5,6), column=c(3,4,5,6))) # m <- matrix(runif(100), nrow=10) # tm = matrix(c(0,1,1,0), nrow=2) # bn <- list(bicluster_net(adjacency_matrix=tm, threshold=.5)) # ensemble_biclusters(bn, b, m)
The function calculates how often features or samples occur across all calculated louvain communities
feature_louvain_overlap(overlap_tables, mat)
feature_louvain_overlap(overlap_tables, mat)
overlap_tables |
List of tables as returned by |
mat |
The data matrix used as input for the biclustering algorithms. |
List of tables as returned by attr_overlap
,
extended by a column showing how often elements occur across all tables.
# a = data.frame(type=c("row", "row", "row", "column", "column", "column"), # ID=c(1,2,3,1,2,3), Fraction=c(1,1,1,.5, .5, .5)) # b = data.frame(type=c("row", "row", "row", "column", "column", "column"), # ID=c(3,2,4,1,5,3), Fraction=c(1,1,1,.5, .5, .5)) # inl <- list(a, b) # feature_louvain_overlap(outl, matrix(1:100, nrow=10))
# a = data.frame(type=c("row", "row", "row", "column", "column", "column"), # ID=c(1,2,3,1,2,3), Fraction=c(1,1,1,.5, .5, .5)) # b = data.frame(type=c("row", "row", "row", "column", "column", "column"), # ID=c(3,2,4,1,5,3), Fraction=c(1,1,1,.5, .5, .5)) # inl <- list(a, b) # feature_louvain_overlap(outl, matrix(1:100, nrow=10))
The function computes a co-occurence network, based on the
function full_graph
.
A similarity threshold is calculated using randomized
biclusters (the bicluster size distribution is kept).
feature_network( bics, mat, n_randomizations = 5, n_steps = 100, plot_edge_dist = TRUE, sn_ratio = 1, error_threshold = 0.05, return_plot_data = FALSE, rr = 1, rc = 1, cc = 1, w = 0 )
feature_network( bics, mat, n_randomizations = 5, n_steps = 100, plot_edge_dist = TRUE, sn_ratio = 1, error_threshold = 0.05, return_plot_data = FALSE, rr = 1, rc = 1, cc = 1, w = 0 )
bics |
A list of bicluster objects. |
mat |
The matrix used for biclustering. |
n_randomizations |
The number of randomizations for cut-off estimation. (The mean of all randomizations is used). |
n_steps |
Number of points where the difference between randomizations and the real data is evaluated. |
plot_edge_dist |
Show the plots for threshold estimation. |
sn_ratio |
If |
error_threshold |
If |
return_plot_data |
Please do not use outside of the package. |
rr |
See |
rc |
See |
cc |
See |
w |
See parameter weighting of |
An object of class cooccurrence_net
.
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # feature_network(bics, m)
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # feature_network(bics, m)
validate_bicluster
.Filter a list of bicluster objects, by erasing all biclusters,
that do not fulfill the minimum number of rows and columns.
Utilizes the function validate_bicluster
.
filter_bicluster_size(bics, minRow, minCol)
filter_bicluster_size(bics, minRow, minCol)
bics |
List of bicluster objects. |
minRow |
Minimum number of rows. |
minCol |
Minimum number of columns. |
A filtered list of bicluster objects.
b <- list(bicluster(row=c(1,2), column=c(1,2,3,4)), bicluster(row=c(3,4,5,6), column=c(3,4,5,6))) filter_bicluster_size(b, 3, 3)
b <- list(bicluster(row=c(1,2), column=c(1,2,3,4)), bicluster(row=c(3,4,5,6), column=c(3,4,5,6))) filter_bicluster_size(b, 3, 3)
If the function returns True
, the bicluster is added to the output
list of biclusters.
Every bicluster is validated, before forwarding to the filter-function.
filter_biclusters(bics, mat, filterfun, ...)
filter_biclusters(bics, mat, filterfun, ...)
bics |
A list of valid bicluster objects. |
mat |
Original matrix, that was used for biclustering. |
filterfun |
A function to filter biclusters. Only if the function
returns |
... |
Other parameters forwarded to the |
A filtered list ob bicluster objects with length(returned_list)<=length(bics).
# m <- matrix(runif(100), nrow=10) b <- list(bicluster(row=c(3,4), column=c(3,4)), bicluster(row=c(3,4,5,6), column=c(3,4,5,6)), bicluster(row=c(3,4,5,6), column=c(3,6))) # filter_biclusters(b, m, function(x) sum(x) < 0)
# m <- matrix(runif(100), nrow=10) b <- list(bicluster(row=c(3,4), column=c(3,4)), bicluster(row=c(3,4,5,6), column=c(3,4,5,6)), bicluster(row=c(3,4,5,6), column=c(3,6))) # filter_biclusters(b, m, function(x) sum(x) < 0)
All values below the threshold will be replaces by 0.
filter_matrix(mat, threshold = 1)
filter_matrix(mat, threshold = 1)
mat |
A Numeric matrix. |
threshold |
All values below will be replaces by 0. |
A filtered numeric matrix.
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # filter_matrix(m, threshold=1)
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # filter_matrix(m, threshold=1)
validate_bicluster
).Remove all biclusters from a list, that are identical
or perfect subsets from each other.
Additionally all invalid biclusters are
removed (See validate_bicluster
).
filter_subsets(bics)
filter_subsets(bics)
bics |
A list of bicluster objects |
A list of bicluster objects, where perfects subsets or identical biclusters are deleted.
filter_subsets(list(bicluster(row=c(1,2,3,4), column=c(1,2,3,4)), bicluster(row=c(1,2,3,4), column=c(1,2,3,4))))
filter_subsets(list(bicluster(row=c(1,2,3,4), column=c(1,2,3,4)), bicluster(row=c(1,2,3,4), column=c(1,2,3,4))))
The function computes a adjacency matrix for rows and columns of biclusters. The matrix values show, how often two rows or two columns or a row and a column occur together in biclusters. In the resulting adjacency matrix, rows are listed first, followed by columns. They have the same order as the the rows and columns of the input matrix.
full_graph( bics, m, rr_weight = 1L, rc_weight = 1L, cc_weight = 1L, weighting = 0L )
full_graph( bics, m, rr_weight = 1L, rc_weight = 1L, cc_weight = 1L, weighting = 0L )
bics |
A list of biclusters. |
m |
The matrix, that was used to calculated the biclusters. |
rr_weight |
Weight row-row interactions. |
rc_weight |
Weight row-col interactions. |
cc_weight |
Weight col-col interactions. |
weighting |
Weight interactions by bicluster size. 0 - no weighting, 1 - multiply by bicluster size, 2 - divide by bicluster size. |
In case the given biclusters have overall more or less columns than rows, the interactions can be weighted to visualize the result properly.
An adjacency matrix.
m <- matrix(seq(1:16), nrow=4) b <- list(bicluster(row=c(1,2,3,4), column=c(1,2,3,4)), bicluster(row=c(3,4,5,6), column=c(3,4,5,6)), bicluster(row=c(3,4,5,6), column=c(3,4,5,6))) # full_graph(b, m)
m <- matrix(seq(1:16), nrow=4) b <- list(bicluster(row=c(1,2,3,4), column=c(1,2,3,4)), bicluster(row=c(3,4,5,6), column=c(3,4,5,6)), bicluster(row=c(3,4,5,6), column=c(3,4,5,6))) # full_graph(b, m)
Return Adjacency matrix from bicluster network
get_adjacency(bic_net)
get_adjacency(bic_net)
bic_net |
An object of class |
Raw unfiltered adjacency matrix.
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # bn <- bicluster_network(bics, m) # get_adjacency(bn)
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # bn <- bicluster_network(bics, m) # get_adjacency(bn)
Return Adjacency matrix from bicluster network
## S4 method for signature 'bicluster_net' get_adjacency(bic_net)
## S4 method for signature 'bicluster_net' get_adjacency(bic_net)
bic_net |
An object of class |
Raw unfiltered adjacency matrix.
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # bn <- bicluster_network(bics, m) # get_adjacency(bn)
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # bn <- bicluster_network(bics, m) # get_adjacency(bn)
Get a unique vector of algorithms from a
list of bicluster
objects.
get_algorithms(bics)
get_algorithms(bics)
bics |
a list of |
A character vector with algorithm names
b <- list(bicluster(row=c(1,2,3,4), column=c(1,2,3,4), algorithm="isa"), bicluster(row=c(3,4,5,6), column=c(3,4,5,6), algorithm="QUBIC"))
b <- list(bicluster(row=c(1,2,3,4), column=c(1,2,3,4), algorithm="isa"), bicluster(row=c(3,4,5,6), column=c(3,4,5,6), algorithm="QUBIC"))
Return algorithms from bicluster network
get_bic_net_algorithms(bic_net)
get_bic_net_algorithms(bic_net)
bic_net |
An object of class |
Algorithm names as characters.
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # bn <- bicluster_network(bics, m) # get_bic_net_algorithms(bn)
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # bn <- bicluster_network(bics, m) # get_bic_net_algorithms(bn)
Return algorithms from bicluster network
## S4 method for signature 'bicluster_net' get_bic_net_algorithms(bic_net)
## S4 method for signature 'bicluster_net' get_bic_net_algorithms(bic_net)
bic_net |
An object of class |
Algorithm names as characters.
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # bn <- bicluster_network(bics, m) # get_bic_net_algorithms(bn)
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # bn <- bicluster_network(bics, m) # get_bic_net_algorithms(bn)
Converts biclusters output of different algorithms/packages in to lists
of bicluster objects. Many algoritms can be directly executed using
the run_...
methods from this package.
This directly returns the converted results. Not all algorithms are shipped
with this package, like Bi-Force, which is running in Java as a standalone
tool or BicARE,
which required an full import using library(BicARE)
in order to run.
But their results can be converted using this function.
get_biclusters(bics, mat, method, transposed = FALSE, filterfun = NULL, ...)
get_biclusters(bics, mat, method, transposed = FALSE, filterfun = NULL, ...)
bics |
A resulting object from a biclustering algorithm (extracted biclusters for fabia) or filename for stored biclustering results. |
mat |
Original matrix, that was used for biclustering. |
method |
Used biclustering package. One of "biclust" (can be further specified as "biclust-bimax", "biclust-cc", "biclust-plaid", "biclust-quest", "biclust-qubic", "biclust-spectral", "biclust-xmotifs", "biclust-unibic"), "BicARE", "isa", "fabia", "biforce", "biclustpy", "qubic2" or "akmbiclust". |
transposed |
Indicate, whether a transposed version of the matrix is
used for biclustering. The |
filterfun |
A function to filter biclusters.
Only if the function returns |
... |
Other parameters forwarded to the |
A list of bicluster
objects, which are
valid (See validate_bicluster
).
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # res <- isa2::isa(m) # get_biclusters(res, m, "isa")
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # res <- isa2::isa(m) # get_biclusters(res, m, "isa")
Extracts the louvain communities from a bicluster_net
or cooccurrence_net
object using the louvain modularity
optimization from the igraph
package (cluster_louvain
).
get_louvain_communities(bic_net, min_size = 2, bics = NULL)
get_louvain_communities(bic_net, min_size = 2, bics = NULL)
bic_net |
A |
min_size |
Minimum size of a louvain community to be returned (minimum value is 2). |
bics |
Optional. Is only use for the class |
A list of bicluster_net
or cooccurrence_net
objects.
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # bn <- bicluster_network(bics, m) # get_louvain_communities(bn)
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # bn <- bicluster_network(bics, m) # get_louvain_communities(bn)
Extracts the louvain communities from a bicluster_net
object using the louvain modularity optimization from the igraph
package (cluster_louvain
).
## S4 method for signature 'bicluster_net' get_louvain_communities(bic_net, min_size = 2, bics = NULL)
## S4 method for signature 'bicluster_net' get_louvain_communities(bic_net, min_size = 2, bics = NULL)
bic_net |
A |
min_size |
Minimum size of a louvain community to be returned. |
bics |
Optional. The respective list of biclusters to identify, from which algorithms a community originates. |
A list of bicluster_net
objects.
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # bn <- bicluster_network(bics, m) # get_louvain_communities(bn)
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # bn <- bicluster_network(bics, m) # get_louvain_communities(bn)
Extracts the louvain communities from a cooccurrence_net
object using the louvain modularity optimization from the igraph
package (cluster_louvain
).
## S4 method for signature 'cooccurrence_net' get_louvain_communities(bic_net, min_size = 2, bics = NULL)
## S4 method for signature 'cooccurrence_net' get_louvain_communities(bic_net, min_size = 2, bics = NULL)
bic_net |
A |
min_size |
Minimum size of a louvain community to be returned (minimum value is 2). |
bics |
Not used. |
A list of cooccurrence_net
objects.
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # fn <- feature_network(bics, m) # get_louvain_communities(fn)
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # fn <- feature_network(bics, m) # get_louvain_communities(fn)
Extract a list of bicluster objects from an akmbiclust biclustering object.
getAkmbiclustClusters(bics, mat, transposed = FALSE, filterfun = NULL, ...)
getAkmbiclustClusters(bics, mat, transposed = FALSE, filterfun = NULL, ...)
bics |
A result object from akmbiclust. |
mat |
Original matrix, that was used for biclustering. |
transposed |
|
filterfun |
A function to filter biclusters.
Only if the function returns |
... |
Other parameters forwarded to the |
A list of bicluster
objects, which have to be
valid (See validate_bicluster
.
# Function called in m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # Not run: run_akmbiclust(m, k=10)
# Function called in m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # Not run: run_akmbiclust(m, k=10)
Get all biclusters from a Bi-Force output file.
getallBFClusters(filename)
getallBFClusters(filename)
filename |
Name of the Bi-Force output file. |
List of biclusters in the form of getBFCluster
a <- "PathToBiForceOutput.txt" # getallBFClusters(a)
a <- "PathToBiForceOutput.txt" # getallBFClusters(a)
Get a bicluster a Bi-Force output file
getBFCluster(filename, cluster)
getBFCluster(filename, cluster)
filename |
Name of the Bi-Force output file. |
cluster |
Number of the bicluster that should be extracted. |
Bicluster as list with rownames in attribute "row" and colnames in attribute "column".
a <- "PathToBiForceOutput.txt" # getBFCluster(a, cluster=1)
a <- "PathToBiForceOutput.txt" # getBFCluster(a, cluster=1)
Extract a list of bicluster objects from an BicARE biclustering object.
getBicAREbiclusters(bics, mat, transposed = FALSE, filterfun = NULL, ...)
getBicAREbiclusters(bics, mat, transposed = FALSE, filterfun = NULL, ...)
bics |
A BicARE bicluster object. |
mat |
Original matrix, that was used for biclustering. |
transposed |
|
filterfun |
A function to filter biclusters. Only if the
function returns |
... |
Other parameters forwarded to the |
A list of bicluster
objects, which have to
be valid (See validate_bicluster
.
# Note that BicARE packackage is not included in the mosbi package m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # res <- BicARE::FLOC(m) # getBicAREbiclusters(res, m)
# Note that BicARE packackage is not included in the mosbi package m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # res <- BicARE::FLOC(m) # getBicAREbiclusters(res, m)
Extract a list of bicluster objects from a biclust object.
getBiclustClusters( bics, mat, method = "biclust", transposed = FALSE, filterfun = NULL, ... )
getBiclustClusters( bics, mat, method = "biclust", transposed = FALSE, filterfun = NULL, ... )
bics |
A biclust object. |
mat |
Original matrix, that was used for biclustering. |
method |
Name of the used biclustering algorithm. Should be one of the following: "biclust", "biclust-bimax", "biclust-cc", "biclust-plaid", "biclust-quest", "biclust-spectral", "biclust-xmotifs" or "biclust-qubic", "biclust-unibic". |
transposed |
|
filterfun |
A function to filter biclusters. Only if the function
returns |
... |
Other parameters forwarded to the |
A list of bicluster
objects, which have
to be valid (See validate_bicluster
.
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # res <- biclust::biclust(m, method = biclust::BCBimax()) # getBiclustClusters(res, m)
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # res <- biclust::biclust(m, method = biclust::BCBimax()) # getBiclustClusters(res, m)
Extract a list of bicluster objects from an biclustpy output file.
getBiclustpyClusters(bics, mat, transposed = FALSE, filterfun = NULL, ...)
getBiclustpyClusters(bics, mat, transposed = FALSE, filterfun = NULL, ...)
bics |
A biclust object. |
mat |
Original matrix, that was used for biclustering. |
transposed |
|
filterfun |
A function to filter biclusters. Only if the
function returns |
... |
Other parameters forwarded to the |
A list of bicluster
objects, which have
to be valid (See validate_bicluster
.
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # Not run: getBiclustpyClusters("PathToFileOfBiclustpyResults", m)
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # Not run: getBiclustpyClusters("PathToFileOfBiclustpyResults", m)
Extract a list of bicluster objects from an fabia biclustering object.
getFabiaClusters(bics, mat, transposed = FALSE, filterfun = NULL, ...)
getFabiaClusters(bics, mat, transposed = FALSE, filterfun = NULL, ...)
bics |
Extracted fabia biclusters. |
mat |
Original matrix, that was used for biclustering. |
transposed |
|
filterfun |
A function to filter biclusters. Only if the function
returns |
... |
Other parameters forwarded to the |
A list of bicluster
objects, which have to
be valid (See validate_bicluster
.
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # res <- fabia::extractBic(fabia::fabia(m, p=5)) # getFabiaClusters(res, m)
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # res <- fabia::extractBic(fabia::fabia(m, p=5)) # getFabiaClusters(res, m)
Extract a list of bicluster objects from an isa2 biclustering object.
getIsaClusters(bics, mat, transposed = FALSE, filterfun = NULL, ...)
getIsaClusters(bics, mat, transposed = FALSE, filterfun = NULL, ...)
bics |
A biclust object. |
mat |
Original matrix, that was used for biclustering. |
transposed |
|
filterfun |
A function to filter biclusters. Only if the function
returns |
... |
Other parameters forwarded to the |
A list of bicluster
objects, which have
to be valid (See validate_bicluster
.
m <- matrix(seq(1:16), nrow=4) # Function part of: # m <- matrix(rnorm(10000), nrow=100) # Not run: run_isa(m)
m <- matrix(seq(1:16), nrow=4) # Function part of: # m <- matrix(rnorm(10000), nrow=100) # Not run: run_isa(m)
Extract biclusters from a QUBIC2 "*.blocks" file. Row and column names are not added to the bicluster objects.
getQUBIC2biclusters(filename, transposed = FALSE)
getQUBIC2biclusters(filename, transposed = FALSE)
filename |
Path to the QUBIC2 results file. |
transposed |
Set to TRUE, if the biclustering was performed on a tranposed matrix. |
A list of validated bicluster
objects (See validate_bicluster
).
a <- "PathToQUBIC2output.txt" # Not run: getQUBIC2biclusters(a)
a <- "PathToQUBIC2output.txt" # Not run: getQUBIC2biclusters(a)
Check, whether a matrix has row- and colnames.
has_names(m)
has_names(m)
m |
A matrix |
Logical indicating existence of row- and colnames.
has_names(matrix(c(1,2,3,4), nrow=2)) m <- matrix(c(1,2,3,4), nrow=2) rownames(m) <- c("r1", "r2") rownames(m) <- c("c1", "c2") has_names(m)
has_names(matrix(c(1,2,3,4), nrow=2)) m <- matrix(c(1,2,3,4), nrow=2) rownames(m) <- c("r1", "r2") rownames(m) <- c("c1", "c2") has_names(m)
Check if a bicluster is a subset (in rows AND columns) of identical to another bicluster.
is_subset_or_identical(bic1, bic2)
is_subset_or_identical(bic1, bic2)
bic1 |
A bicluster. |
bic2 |
A bicluster. |
1 if bic1 is a subset of bic2, 2 if bic 1 is identical to bic2, 0 else.
is_subset_or_identical(bicluster(row=c(1,2,3,4), column=c(1,2,3,4)), bicluster(row=c(1,2,3,4), column=c(1,2,3,4)))
is_subset_or_identical(bicluster(row=c(1,2,3,4), column=c(1,2,3,4)), bicluster(row=c(1,2,3,4), column=c(1,2,3,4)))
Development of the mice brain lipidome over several weeks.
data(mouse_data)
data(mouse_data)
A data frame with 245 rows and 61 columns
https://www.ebi.ac.uk/metabolights/MTBLS562
# data(mouse_data) # View(mouse_data)
# data(mouse_data) # View(mouse_data)
Computes the how many edges remain in a network if edges with a weight
lower than a certain threshold are removed.
The number of remaining edges between 1 and max(adjm) are calculated.
It is assumend that the matrix is symmetric and therefore the number
of edges divided by two.
Uses the function replace_values
.
network_edge_strength(adjm)
network_edge_strength(adjm)
adjm |
A symmetrix numeric matrix. |
A numeirc matrix of dim(max(adjm), 2)
. The first column
indicated the applied threshold, the second column the remaining edges.
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # fn <- feature_network(bics, m) # network_edge_strength(apply_threshold(fn))
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # fn <- feature_network(bics, m) # network_edge_strength(apply_threshold(fn))
Same as network_edge_strength
, but for (positive)
non-integer matrices.
network_edge_strength_float(adjm, steps = 100L, maximum = 0)
network_edge_strength_float(adjm, steps = 100L, maximum = 0)
adjm |
A symmetrix numeric matrix. |
steps |
Number of steps for which the edge count is evaluated. |
maximum |
Highest value until which the edge weight is evaluated.
If maximum=0, the max value of |
Computes the how many edges remain in a network if edges with a
weight lower than a certain threshold are removed.
The number of remaining edges between 1 and max(adjm) are calculated.
It is assumend that the matrix is symmetric and therefore
the number of edges divided by two.
Uses the function replace_values_float
.
A numeirc matrix of dim(max(adjm), 2)
.
The first column indicated the applied threshold, the second
column the remaining edges.
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # bn <- bicluster_network(bics, m) # network_edge_strength_float(apply_threshold(bn))
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # bn <- bicluster_network(bics, m) # network_edge_strength_float(apply_threshold(bn))
Get the number of biclusters, generated by the Bi-Force algorithm.
NoBFBiclusters(filename)
NoBFBiclusters(filename)
filename |
Name of the Bi-Force output file. |
Number of biclusters.
a <- "PathToBiForceOutput.txt" # NoBFBiclusters(a)
a <- "PathToBiForceOutput.txt" # NoBFBiclusters(a)
When plotting bicluster networks, node sizes adapted to bicluster sizes
can improve visual inspection.
Node sizes are computed using the following formula:
(atan( (x - min(x)) / (max(x) - min(x)) + offset ) * base_size)
.
With x being defined a vector of bicluster sizes defined
by the MARGIN
parameter.
node_size(bics, base_size = 10, offset = 0.2, MARGIN = "column")
node_size(bics, base_size = 10, offset = 0.2, MARGIN = "column")
bics |
A list of |
base_size |
Is multiplied with the atan result for the node size |
offset |
Offset for the atan calculation. Has to be > 0. Smaller values result in higher differences of node sizes. |
MARGIN |
"column", "row" or "both" are taken into account for the size of a bicluster bicluster |
Vector of node sizes as floats.
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # bn <- bicluster_network(bics, m) # nz <- node_size(bics) # plot_algo_network(bn, bics, vertex.size=nz) # plot(bn, vertex.size=node_size(bics, offset=.1, base_size=15))
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # bn <- bicluster_network(bics, m) # nz <- node_size(bics) # plot_algo_network(bn, bics, vertex.size=nz) # plot(bn, vertex.size=node_size(bics, offset=.1, base_size=15))
The function computes a matrix with the same dimensions as the input matrix and fills the matrix elements with the frequence of occurance of the data points in the input list of biclusters.
occurance_matrix(bics, mat)
occurance_matrix(bics, mat)
bics |
A list of |
mat |
The data matrix used for biclustering. |
A numeric matrix with the dimensions of the input matrix. The values represent the frequency of occurance of this point in the list of biclusters.
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # occurance_matrix(bics, m)
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # occurance_matrix(bics, m)
The function uses the occurance_matrix
function and
returns all values higher than the threshold
as a DataFrame.
occurance_table(bics, mat, threshold = 0)
occurance_table(bics, mat, threshold = 0)
bics |
A list of |
mat |
The data matrix used for biclustering. |
threshold |
Only data points higher than this threshold are returned. |
A DataFrame with the frequencies of occurance for values higher
than a threshold
.
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # occurance_table(bics, m, threshold=.1)
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # occurance_table(bics, m, threshold=.1)
The probability is computed using the
forumla .
p_overlap(x, y, k, n)
p_overlap(x, y, k, n)
x |
Overlap. |
y |
Size of sample 1. |
k |
Size of Sample 2. |
n |
Number of all elements sampled from. |
Overlap probability.
p_overlap(10, 20, 30, 100)
p_overlap(10, 20, 30, 100)
Is computed by calculating the overlap probability for each
dimension independently and multiplying them using the
function p_overlap
.
p_overlap_2d(ov_x, ov_y, s1x, s1y, s2x, s2y, mat_x, mat_y)
p_overlap_2d(ov_x, ov_y, s1x, s1y, s2x, s2y, mat_x, mat_y)
ov_x |
Overlap in the first dimension. |
ov_y |
Overlap in the second dimension. |
s1x |
First sample of the first dimension. |
s1y |
First sample of the second dimension. |
s2x |
Second sample of first dimension. |
s2y |
Second sample of the second dimension. |
mat_x |
Number of all elements from the first dimension sampled from. |
mat_y |
Number of all elements from the second dimension sampled from. |
Overlap probability.
p_overlap_2d(10, 10, 20, 20, 30, 30, 100, 100)
p_overlap_2d(10, 10, 20, 20, 30, 30, 100, 100)
Is computed by adding up probabilities for all combinations of
the observed or higher overlaps using the
function p_overlap_2d
.
p_overlap_2d_higher(ov_x, ov_y, s1x, s1y, s2x, s2y, mat_x, mat_y)
p_overlap_2d_higher(ov_x, ov_y, s1x, s1y, s2x, s2y, mat_x, mat_y)
ov_x |
Overlap in the first dimension. |
ov_y |
Overlap in the second dimension. |
s1x |
First sample of the first dimension. |
s1y |
First sample of the second dimension. |
s2x |
Second sample of first dimension. |
s2y |
Second sample of the second dimension. |
mat_x |
Number of all elements from the first dimension sampled from. |
mat_y |
Number of all elements from the second dimension sampled from. |
Overlap probability
p_overlap_2d_higher(10, 10, 20, 20, 30, 30, 100, 100)
p_overlap_2d_higher(10, 10, 20, 20, 30, 30, 100, 100)
Is computed by adding up probabilities for all possible
overlaps equal or higher to the observed one using the
function p_overlap
.
p_overlap_higher(x, y, k, n)
p_overlap_higher(x, y, k, n)
x |
Overlap. |
y |
Size of sample 1. |
k |
Size of Sample 2. |
n |
Number of all elements sampled from. |
Overlap probability.
p_overlap_higher(10, 20, 30, 100)
p_overlap_higher(10, 20, 30, 100)
In the plot each bicluster is colored by the algorithm, that generated it.
plot_algo_network(bic_net, bics, new_layout = TRUE, ...)
plot_algo_network(bic_net, bics, new_layout = TRUE, ...)
bic_net |
A |
bics |
The corresponding list of biclusters from |
new_layout |
If |
... |
Plot parameters forwarded
to |
If new_layout
, a new network layout is returned
that can be used for other plots.
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # bn <- bicluster_network(bics, m) # plot_algo_network(bn, bics)
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # bn <- bicluster_network(bics, m) # plot_algo_network(bn, bics)
Plot a bicluster network with piecharts as nodes.
plot_piechart_bicluster_network( bic_net, bics, class_vector, colors, named = TRUE, MARGIN = "column", new_layout = TRUE, ... )
plot_piechart_bicluster_network( bic_net, bics, class_vector, colors, named = TRUE, MARGIN = "column", new_layout = TRUE, ... )
bic_net |
A |
bics |
The corresponding list of biclusters from |
class_vector |
A (named) vector with class affinities. Every occuring element in the biclustes must have a non NA value in this list. |
colors |
Colors used for the classes. Must be a vector with colors in the order of sort(unique(class_vector)). |
named |
Indicates if |
MARGIN |
Must be "row" or "column". Indicates which dimension of the bicluster should be used for coloring. |
new_layout |
If |
... |
Additional parameters forwarded
to |
If new_layout
, a new network layout is returned that can
be used for other plots.
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # bn <- bicluster_network(bics, m) # groups <- ifelse(runif(100)< 0.5, "group1", "group2") # cols <- c("group1"="blue", "group2"="grey") # plot_piechart_bicluster_network(bn, bics, groups, cols, named=FALSE)
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # bn <- bicluster_network(bics, m) # groups <- ifelse(runif(100)< 0.5, "group1", "group2") # cols <- c("group1"="blue", "group2"="grey") # plot_piechart_bicluster_network(bn, bics, groups, cols, named=FALSE)
Converts the object into a graph and uses its plot function.
## S4 method for signature 'bicluster_net,missing' plot(x, y, ...)
## S4 method for signature 'bicluster_net,missing' plot(x, y, ...)
x |
An object of class |
y |
Not used. |
... |
Plot parameters forwarded
to |
An graph plot.
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # bn <- bicluster_network(bics, m) # plot(bn)
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # bn <- bicluster_network(bics, m) # plot(bn)
Converts the object into a graph and uses its plot function.
## S4 method for signature 'cooccurrence_net,missing' plot(x, y, ...)
## S4 method for signature 'cooccurrence_net,missing' plot(x, y, ...)
x |
An object of class |
y |
Not used. |
... |
Plot parameters forwarded
to |
An graph plot.
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # fn <- feature_network(bics, m) # plot(fn)
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # fn <- feature_network(bics, m) # plot(fn)
Randomize a matrix bu shuffling all rows and columns.
randomize_matrix(m)
randomize_matrix(m)
m |
A matrix. |
A randomized version of the input matrix.
m <- matrix(c(1,2,3,4), nrow=2) randomize_matrix(m)
m <- matrix(c(1,2,3,4), nrow=2) randomize_matrix(m)
This function replaces all elements of an integer matrix, which are under a certain threshold (<) with zero.
replace_threshold(m, threshold)
replace_threshold(m, threshold)
m |
A numeric matrix. |
threshold |
A numeric threshold under which all elements in the matrix are replaced by zero. |
An integer matrix.
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # replace_threshold(m, 1)
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # replace_threshold(m, 1)
Replace values in an integer matrix, that are lower than a certain threshold.
replace_values(mat, threshold, replace_higher = TRUE)
replace_values(mat, threshold, replace_higher = TRUE)
mat |
An integer matrix |
threshold |
All values in the matrix lower than this values are replaced by 0. |
replace_higher |
If set to true, all values >= |
An integer matrix with (partially) replaced values.
replace_values(matrix(seq(1, 16), nrow=4), threshold=4)
replace_values(matrix(seq(1, 16), nrow=4), threshold=4)
Same as replace_values
, but for (positive) non-integer
matrices.
replace_values_float(mat, threshold, replace_higher = TRUE)
replace_values_float(mat, threshold, replace_higher = TRUE)
mat |
A numeric matrix |
threshold |
All values in the matrix lower than this values are replaced by 0. |
replace_higher |
If set to true, all values >= |
Replace values in a numeric matrix, that are lower than a certain threshold.
A numeric matrix with (partially) replaced values.
replace_values(matrix(rnorm(100), nrow=10), threshold=1)
replace_values(matrix(rnorm(100), nrow=10), threshold=1)
Can be used for e.g. histograms.
rowhistogram(bic)
rowhistogram(bic)
bic |
A list of bicluster objects. |
A vector with the lenghts of the rows in every bicluster object.
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # rowhistogram(bics)
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # rowhistogram(bics)
The function executes the akmbiclust biclustering algorithm, returning a list of biclusters converted into bicluster objects compatible with this package. If the algorithm fails to run, an empty list is returned.
run_akmbiclust(data_matrix, minRow = 2, minCol = 2, ...)
run_akmbiclust(data_matrix, minRow = 2, minCol = 2, ...)
data_matrix |
A numeric matrix. |
minRow |
Same parameters as in filter_bicluster_size. |
minCol |
Same parameters as in filter_bicluster_size. |
... |
Other parameters forwarded to the akmbiclust function. |
a list of bicluster objects.
m <- matrix(seq(1:16), nrow=4) # set.seed(10) # m <- matrix(rnorm(10000), nrow=100) # Not run: run_akmbiclust(m, k=10)
m <- matrix(seq(1:16), nrow=4) # set.seed(10) # m <- matrix(rnorm(10000), nrow=100) # Not run: run_akmbiclust(m, k=10)
The function executes the BCBimax biclustering algorithm, returning a list of biclusters converted into bicluster objects compatible with this package. If the algorithm fails to run, an empty list is returned.
run_bimax(data_matrix, minRow = 2, minCol = 2, ...)
run_bimax(data_matrix, minRow = 2, minCol = 2, ...)
data_matrix |
A numeric matrix. |
minRow |
Same parameters as in filter_bicluster_size. |
minCol |
Same parameters as in filter_bicluster_size. |
... |
Other parameters forwarded to the BCBimax function. |
a list of bicluster objects.
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # run_bimax(m)
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # run_bimax(m)
The function executes the BCCC biclustering algorithm, returning a list of biclusters converted into bicluster objects compatible with this package. If the algorithm fails to run, an empty list is returned.
run_cc(data_matrix, minRow = 2, minCol = 2, ...)
run_cc(data_matrix, minRow = 2, minCol = 2, ...)
data_matrix |
A numeric matrix. |
minRow |
Same parameters as in filter_bicluster_size. |
minCol |
Same parameters as in filter_bicluster_size. |
... |
Other parameters forwarded to the BCCC function. |
a list of bicluster objects.
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # run_cc(m)
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # run_cc(m)
The function executes the fabia biclustering algorithm, returning a list of biclusters converted into bicluster objects compatible with this package. If the algorithm fails to run, an empty list is returned.
run_fabia( data_matrix, minRow = 2, minCol = 2, thresZ = 0.5, thresL = NULL, ... )
run_fabia( data_matrix, minRow = 2, minCol = 2, thresZ = 0.5, thresL = NULL, ... )
data_matrix |
A numeric matrix. |
minRow |
Same parameters as in filter_bicluster_size. |
minCol |
Same parameters as in filter_bicluster_size. |
thresZ |
See parameter from the extractBic function. |
thresL |
See parameter from the extractBic function. |
... |
Other parameters forwarded to the fabia function. |
a list of bicluster objects.
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(1000), nrow=10) # run_fabia(m, p=5)
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(1000), nrow=10) # run_fabia(m, p=5)
The function executes the isa biclustering algorithm, returning a list of biclusters converted into bicluster objects compatible with this package. If the algorithm fails to run, an empty list is returned.
run_isa(data_matrix, minRow = 2, minCol = 2, ...)
run_isa(data_matrix, minRow = 2, minCol = 2, ...)
data_matrix |
A numeric matrix. |
minRow |
Same parameters as in filter_bicluster_size. |
minCol |
Same parameters as in filter_bicluster_size. |
... |
Other parameters forwarded to the isa function. |
a list of bicluster objects.
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # Not run: run_isa(m)
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # Not run: run_isa(m)
The function executes the BCPlaid biclustering algorithm, returning a list of biclusters converted into bicluster objects compatible with this package. If the algorithm fails to run, an empty list is returned.
run_plaid(data_matrix, minRow = 2, minCol = 2, ...)
run_plaid(data_matrix, minRow = 2, minCol = 2, ...)
data_matrix |
A numeric matrix. |
minRow |
Same parameters as in filter_bicluster_size. |
minCol |
Same parameters as in filter_bicluster_size. |
... |
Other parameters forwarded to the BCPlaid function. |
a list of bicluster objects.
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # run_plaid(m)
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # run_plaid(m)
The function executes the BCQU biclustering algorithm, returning a list of biclusters converted into bicluster objects compatible with this package. If the algorithm fails to run, an empty list is returned.
run_qubic(data_matrix, minRow = 2, minCol = 2, ...)
run_qubic(data_matrix, minRow = 2, minCol = 2, ...)
data_matrix |
A numeric matrix. |
minRow |
Same parameters as in filter_bicluster_size. |
minCol |
Same parameters as in filter_bicluster_size. |
... |
Other parameters forwarded to the BCQU function. |
a list of bicluster objects.
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # run_qubic(m)
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # run_qubic(m)
The function executes the BCQuest biclustering algorithm, returning a list of biclusters converted into bicluster objects compatible with this package. If the algorithm fails to run, an empty list is returned.
run_quest(data_matrix, minRow = 2, minCol = 2, ...)
run_quest(data_matrix, minRow = 2, minCol = 2, ...)
data_matrix |
A numeric matrix. |
minRow |
Same parameters as in filter_bicluster_size. |
minCol |
Same parameters as in filter_bicluster_size. |
... |
Other parameters forwarded to the BCQuest function. |
a list of bicluster objects.
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # run_quest(m)
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # run_quest(m)
The function executes the BCSpectral biclustering algorithm, returning a list of biclusters converted into bicluster objects compatible with this package. If the algorithm fails to run, an empty list is returned.
run_spectral(data_matrix, minRow = 2, minCol = 2, ...)
run_spectral(data_matrix, minRow = 2, minCol = 2, ...)
data_matrix |
A numeric matrix. |
minRow |
Same parameters as in filter_bicluster_size. |
minCol |
Same parameters as in filter_bicluster_size. |
... |
Other parameters forwarded to the BCSpectral function. |
a list of bicluster objects.
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # run_spectral(m)
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # run_spectral(m)
The function executes the runibic::BCUnibic biclustering algorithm, returning a list of biclusters converted into bicluster objects compatible with this package. If the algorithm fails to run, an empty list is returned.
run_unibic(data_matrix, minRow = 2, minCol = 2, ...)
run_unibic(data_matrix, minRow = 2, minCol = 2, ...)
data_matrix |
A numeric matrix. |
minRow |
Same parameters as in filter_bicluster_size. |
minCol |
Same parameters as in filter_bicluster_size. |
... |
Other parameters forwarded to the runibic::BCUnibic function. function. |
a list of bicluster objects.
Function as a string, which can be executed.
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # run_unibic(m, nbic=10)
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # run_unibic(m, nbic=10)
The function executes the BCXmotifs biclustering algorithm, returning a list of biclusters converted into bicluster objects compatible with this package. If the algorithm fails to run, an empty list is returned.
run_xmotifs(data_matrix, minRow = 2, minCol = 2, ...)
run_xmotifs(data_matrix, minRow = 2, minCol = 2, ...)
data_matrix |
A numeric matrix. |
minRow |
Same parameters as in filter_bicluster_size. |
minCol |
Same parameters as in filter_bicluster_size. |
... |
Other parameters forwarded to the BCXmotifs function. |
a list of bicluster objects.
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # run_xmotifs(m)
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # run_xmotifs(m)
The function generates a list of biclusters given an input list of biclusters, where each bicluster has the same number or rows and columns, but with sampled entries from a uniform distribution of all rows and columns is the matrix.
sample_biclusters(bics, mat)
sample_biclusters(bics, mat)
bics |
A list of validated bicluster objects. |
mat |
The numeric matrix, that was used to generate the biclusters. |
A list of bicluster objects.
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # sample_biclusters(bics, m)
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # sample_biclusters(bics, m)
The function returns an adapted bicluster list based on
a bicluster_net
object.
This might be necessary e.g. after get_louvain_communities
was used a community consists only of a subset of the biclusters.
select_biclusters_from_bicluster_network(bic_net, bics)
select_biclusters_from_bicluster_network(bic_net, bics)
bic_net |
|
bics |
A list of |
A subsetted list of bicluster
objects
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # bn <- bicluster_network(bics, m) # lc <- get_louvain_communities(bn) # select_biclusters_from_bicluster_network(lc[[1]], bics)
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # bn <- bicluster_network(bics, m) # lc <- get_louvain_communities(bn) # select_biclusters_from_bicluster_network(lc[[1]], bics)
The function returns an adapted bicluster list based on
a bicluster_net
object.
This might be necessary e.g. after get_louvain_communities
was used and a community consists only of a subset of the biclusters.
## S4 method for signature 'bicluster_net,list' select_biclusters_from_bicluster_network(bic_net, bics)
## S4 method for signature 'bicluster_net,list' select_biclusters_from_bicluster_network(bic_net, bics)
bic_net |
|
bics |
A list of |
A subsetted list of bicluster
objects.
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # bn <- bicluster_network(bics, m) # lc <- get_louvain_communities(bn) # select_biclusters_from_bicluster_network(lc[[1]], bics)
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # bn <- bicluster_network(bics, m) # lc <- get_louvain_communities(bn) # select_biclusters_from_bicluster_network(lc[[1]], bics)
Add row-/colnames to a bicluster object.
set_bicluster_names(bic, m)
set_bicluster_names(bic, m)
bic |
A bicluster object. |
m |
The matrix, that was used for the biclustering. (Works only if matrix has row-/colnames.) |
The updated bicluster object.
m <- matrix(c(1,2,3,4), nrow=2) rownames(m) <- c("r1", "r2") rownames(m) <- c("c1", "c2") set_bicluster_names(bicluster(row=c(1,2), column=c(1,2)), m)
m <- matrix(c(1,2,3,4), nrow=2) rownames(m) <- c("r1", "r2") rownames(m) <- c("c1", "c2") set_bicluster_names(bicluster(row=c(1,2), column=c(1,2)), m)
Add row-/colnames to a bicluster object.
## S4 method for signature 'bicluster,matrix' set_bicluster_names(bic, m)
## S4 method for signature 'bicluster,matrix' set_bicluster_names(bic, m)
bic |
A bicluster object. |
m |
The matrix, that was used for the biclustering. (Works only if matrix has row-/colnames.) |
The updated bicluster object.
#' @examples m <- matrix(c(1,2,3,4), nrow=2) rownames(m) <- c("r1", "r2") rownames(m) <- c("c1", "c2") set_bicluster_names(bicluster(row=c(1,2), column=c(1,2)), m)
This function computes a similarity matrix between biclusters using different similarity metrics.
similarity_matrix( bics, MARGIN = "both", metric = 1L, prob_scale = FALSE, mat_row = 0L, mat_col = 0L, prl = FALSE )
similarity_matrix( bics, MARGIN = "both", metric = 1L, prob_scale = FALSE, mat_row = 0L, mat_col = 0L, prl = FALSE )
bics |
A list of bicluster objects. |
MARGIN |
Choose if the distance is computed over |
metric |
Integer indicating which metric is used. 1: Bray-Curtis similarity (default), 2: Jaccard index, 3: overlap coefficient, 4: Fowlkes–Mallows index. |
prob_scale |
Scale similarity by the probability of an
overlap equal of higher to the observed one. The scaling is
done by multiplying the similarity
with |
mat_row |
If |
mat_col |
If |
prl |
Compute the similarity matrix using multiple
cores (works only for |
A numeric matrix of the similarities between all given biclusters.
b <- list(bicluster(row=c(1,2,3,4), column=c(1,2,3,4)), bicluster(row=c(3,4,5,6), column=c(3,4,5,6))) similarity_matrix(b)
b <- list(bicluster(row=c(1,2,3,4), column=c(1,2,3,4)), bicluster(row=c(3,4,5,6), column=c(3,4,5,6))) similarity_matrix(b)
Transpose a bicluster. Row and column slots will be changed.
transpose_bicluster(bic)
transpose_bicluster(bic)
bic |
A bicluster object. |
A transposed bicluster object,
transpose_bicluster(bicluster(row=c(3,4,5,6), column=c(3,4,5,6)))
transpose_bicluster(bicluster(row=c(3,4,5,6), column=c(3,4,5,6)))
Indicates, whether a bicluster is valid. That means it needs at least one row and one column.
validate_bicluster(bic, minRow = 1L, minCol = 1L)
validate_bicluster(bic, minRow = 1L, minCol = 1L)
bic |
A bicluster object |
minRow |
Minimum number of required rows (Min=1). |
minCol |
Minimum number of required columns (Min=1). |
Logical indicating a valid bicluster object.
validate_bicluster(bicluster(row=c(3,4,5,6), column=c(3,4,5,6)))
validate_bicluster(bicluster(row=c(3,4,5,6), column=c(3,4,5,6)))
Save and adjacency matrix as returned by full_graph
or
1 - distance_matrix
as a GraphML file.
write_graphml(m, filename, cols)
write_graphml(m, filename, cols)
m |
A symmetric numeric matrix (Adjacency matrix). Rownames are considered as node names. |
filename |
Name of the resulting GraphML file (should end with ".gml"). |
cols |
Node colors. |
0 if successful.
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # bn <- bicluster_network(bics, m) # write_graphml(apply_threshold(bn), "testfile.txt")
m <- matrix(seq(1:16), nrow=4) # m <- matrix(rnorm(10000), nrow=100) # bics <- c(run_fabia(m), run_isa(m), run_plaid(m)) # bn <- bicluster_network(bics, m) # write_graphml(apply_threshold(bn), "testfile.txt")
Write an R matrix to a file (In a Bi-Force or QUBIC2 readable format).
write_matrix(m, filename, qubic2_format = FALSE)
write_matrix(m, filename, qubic2_format = FALSE)
m |
A Numeric matrix. |
filename |
Name of the output file. |
qubic2_format |
Write the matrix in a format QUBIC2 is able to read. This means adding a row- and column names to the file. |
0 if file was written successfully.
write_matrix(matrix(c(1,2,3,4), nrow=2), "testfile.txt")
write_matrix(matrix(c(1,2,3,4), nrow=2), "testfile.txt")
Make a vector of R indices compatible with c++ by substracting every element by one.
zero_subsetting(v)
zero_subsetting(v)
v |
A numeric vector. |
A numeric vector with every element decremented by one.
zero_subsetting(c(1,2,3,4,5))
zero_subsetting(c(1,2,3,4,5))