birewire.analysis.bipartite: Analysis of Jaccard similarity trends across switching steps.

Description Usage Arguments Details Value Author(s) References Examples

View source: R/BiRewire.R

Description

This function performs a sequence of max.iter switching steps on the input bipartite graph g and compute the Jaccard similarity between g (the initial network) and its rewired version each step switching steps. This procedure is pefromed n.networks times and a simple explorative plot, with mean and CI, is visualized if display is set to true.

Usage

1
2
birewire.analysis.bipartite(incidence, step=10, max.iter="n",accuracy=0.00005,
verbose=TRUE,MAXITER_MUL=10,exact=FALSE,n.networks=50,display=TRUE)

Arguments

incidence

Incidence matrix of the initial bipartite graph g (can be extracted from an igraph bipartite graph using the get.incidence function). Since 3.6.0 this matrix can contain also NAs and the position of such entries will be preserved by the SA;

step

10 (default): the interval (in terms of switching steps) at which the Jaccard index between g and the its current rewired version is computed;

max.iter

"n" (default) the number of switching steps to be performed (or if exact==TRUE the number of successful switching steps). If equal to "n" then this number is considered equal to the analytically derived lower bound presented in Gobbi et al. (see References): N={e}/{2(1-d)} \ln{((e-de)/δ)} if exact is FALSE, N={e(1-d)}/{2} \ln{((e-de)/δ)} otherwise , where e is the number of edges of g and d its edge density . This bound is much lower than the empirical one proposed in Milo et al. 2003 (see References);

accuracy

0.00005 (default) is the desired level of accuracy reflecting the average distance between the Jaccard index at the N-th step and its analytically derived fixed point in terms of fracion of common edges;

verbose

TRUE (default). When TRUE a progression bar is printed during computation;

MAXITER_MUL

10 (default). If exact==TRUE in order to prevent a possible infinite loop the program stops anyway after MAXITER_MUL*max.iter iterations;

exact

FALSE (default). If TRUE the program performs max.iter swithcing steps, otherwise the program will count also the not-performed swithcing steps;

n.networks

50 (default), the number of independent rewiring process starting from the same inital graph from which the mean value and the CI is computed.

display

TRUE (default). If TRUE two explorative plots are displayed summarizing the trend of the Jaccard index in terms of mean and confidence interval.

Details

This function performs max.iter switching steps (see references). In particular, at each step two edges are randomly selected from the current version of g. Let these two edges be (a,b) and (c,d) (where a and c belong to the first class of nodes whereas b and d belong to the second one), with a\not=c and b\not=d.
If the (a,d) and (c,b) edges are not already present in the current current version of g then (a,d) and (c,b) replace (a,b) and (c,d).

At each step number of switching steps the function computes the Jaccard index between the original graph g and its current version.
This procedure is perfomed n.networks times and if display is set to TRUE, two explorative plots showing the mean value of the Jaccad Index over the SS and its CI are displayed.

Value

A list containing a data.frame data collecting all the Jacard index computed (each row is a run of the SA), and the analytically derived lower bound N of switching steps to be performed by the switching algorithm in order to provide the revired version of g with the maximal level of achievable randomness (in terms of dissimilarity from the initial g).

Author(s)

Andrea Gobbi
Maintainer: Andrea Gobbi <gobbi.andrea@mail.com>
Special thanks to:
Davide Albanese

References

Gobbi, A. and Iorio, F. and Dawson, K. J. and Wedge, D. C. and Tamborero, D. and Alexandrov, L. B. and Lopez-Bigas, N. and Garnett, M. J. and Jurman, G. and Saez-Rodriguez, J. (2014) Fast randomization of large genomic datasets while preserving alteration counts Bioinformatics 2014 30 (17): i617-i623 doi: 10.1093/bioinformatics/btu474.

Iorio, F. and and Bernardo-Faura, M. and Gobbi, A. and Cokelaer, T.and Jurman, G.and Saez-Rodriguez, J. (2016) Efficient randomization of biologicalnetworks while preserving functionalcharacterization of individual nodes Bioinformatics 2016 1 (17):542 doi: 10.1186/s12859-016-1402-1.

Jaccard, P. (1901), Étude comparative de la distribution florale dans une portion des Alpes et des Jura, Bulletin de la Société Vaudoise des Sciences Naturelles 37: 547–579.

R. Milo, N. Kashtan, S. Itzkovitz, M. E. J. Newman, U. Alon (2003), On the uniform generation of random graphs with prescribed degree sequences, eprint arXiv:cond-mat/0312028

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
library(BiRewire)
g <-graph.bipartite(rep(0:1,length=10), c(1:10))

##get the incidence matrix of g
 m<-as.matrix(get.incidence(graph=g))

## set parameters
step=1
max=100*length(E(g))

## perform two different analysis using two different maximal number of switching steps
scores<-birewire.analysis.bipartite(m,step,max,n.networks=10)
scores2<-birewire.analysis.bipartite(m,step,"n",n.networks=10)

BiRewire documentation built on Nov. 8, 2020, 8:09 p.m.