Title: | Time-Ordered and Time-Aggregated Network Analyses |
---|---|
Description: | Approaches for incorporating time into network analysis. Methods include: construction of time-ordered networks (temporal graphs); shortest-time and shortest-path-length analyses; resource spread calculations; data resampling and rarefaction for null model construction; reduction to time-aggregated networks with variable window sizes; application of common descriptive statistics to these networks; vector clock latencies; and plotting functionalities. The package supports <doi:10.1371/journal.pone.0020298>. |
Authors: | Benjamin Wong Blonder [aut, cre] |
Maintainer: | Benjamin Wong Blonder <[email protected]> |
License: | GPL-3 |
Version: | 1.0.1 |
Built: | 2024-11-14 04:59:10 UTC |
Source: | https://github.com/bblonder/timeordered |
From a recent study of information flow in ant colonies. In this study, ants were uniquely marked with paint and identified by a four letter code - e.g. WGWB denotes an ant with a red head, green thorax, white left gaster, and blue right gaster. Body positions with missing paint marks are denoted with underscores.
In-nest activity was recorded with a high definition video camera. The complete set of pairwise interactions between all individuals at all times was obtained by several undergraduates repeatedly watching each video. Interactions were defined as any touch between one ant's antenna and any body part of another ant.
The dataset contains four columns: VertexFrom, VertexTo, TimeStart, and TimeStop. Each row is a unique interaction between two ants. Each interaction is directed, indicating that the VertexFrom ant has initiated a contact with the VertexTo ant. TimeStart and TimeStop characterize when the interaction began and finished. In this demo version of the data set, TimeStop = TimeStart + 1. Times are recorded in seconds.
ants
ants
A data frame containing 1911 observations over 24 minutes.
Blonder & Dornhaus (2011), Supplementary Information, Colony 1-1.
Blonder & Dornhaus, Time-ordered networks reveal limitations to information flow in ant colonies. PLoS One (2011), in press.
-
applynetworkfunction(slices, fun)
applynetworkfunction(slices, fun)
slices |
A list of time-aggregated networks, of class igraph |
fun |
The function to be applied; takes a single argument |
A list whose entries represent the function's value for each network
Benjamin Blonder [email protected].
data(ants) allindivs <- c(union(as.character(ants$VertexFrom), as.character(ants$VertexTo)), "NULL1", "NULL2") g <- generatetonetwork(ants, allindivs) td100 <- generatetimedeltas(0,1500,100) ns100 <- generatenetworkslices(g, td100) md100 <- applynetworkfunction(ns100, diameter) tl100 <- generatetimelags(0,1500,100) nl100 <- generatenetworkslices(g, tl100) ml100 <- applynetworkfunction(nl100, function(x){diameter(x)}) par(mfrow=c(1,2)) plot(midpoints(td100),unlist(md100),type="l",xlab="Time (window size = 100)",ylab="Diameter") plot(maxpoints(tl100),unlist(ml100),type="l",xlab="Aggregation time",ylab="Diameter")
data(ants) allindivs <- c(union(as.character(ants$VertexFrom), as.character(ants$VertexTo)), "NULL1", "NULL2") g <- generatetonetwork(ants, allindivs) td100 <- generatetimedeltas(0,1500,100) ns100 <- generatenetworkslices(g, td100) md100 <- applynetworkfunction(ns100, diameter) tl100 <- generatetimelags(0,1500,100) nl100 <- generatenetworkslices(g, tl100) ml100 <- applynetworkfunction(nl100, function(x){diameter(x)}) par(mfrow=c(1,2)) plot(midpoints(td100),unlist(md100),type="l",xlab="Time (window size = 100)",ylab="Diameter") plot(maxpoints(tl100),unlist(ml100),type="l",xlab="Aggregation time",ylab="Diameter")
Vector clock latencies describe the minimum time delay between one individual broadcasting a signal and another individual receiving it, at a given time, through any causally permitted path in the time-ordered network. Smaller values indicate individuals that are connected by shorter causally-permitted paths at a given time.
generatelatencies(raw, allindivs)
generatelatencies(raw, allindivs)
raw |
An event list, consisting of a data frame with four columns: VertexFrom, VertexTo, TimeStart, and TimeStop. Each row in this data frame represents a single directed interaction event between VertexFrom and VertexTo beginning at TimeStart and ending at TimeStop. Assumes that no event begins at a time less than zero. |
allindivs |
A list of all possible vertices including ones not observed interacting during the range of time reported in |
A n x n x m array, where n is the number of vertices and m is the maximum start time in the raw event list. The [i,j,k] entry of the array describes the latency from i to j at time k. NA is returned if there is not causally permitted path between i and j by time k.
Return value can require large memory allocation depending on the data set. Ensure that data contains no times < 0 before running.
Benjamin Blonder [email protected].
Kossinets et al. The structure of information pathways in a social communication network. KDD '08: Proceeding of the 14th ACM SIGKDD international conference on Knowledge discovery and data mining (2008)
data(ants) allindivs <- c(union(as.character(ants$VertexFrom), as.character(ants$VertexTo)), "NULL1", "NULL2") l <- generatelatencies(ants, allindivs) image(l[,,1000],axes=FALSE,frame=TRUE,col=rainbow(100)) axis(1, at = (1:ncol(l))/ncol(l), labels=colnames(l),tick=FALSE,las=2,cex.axis=0.2) axis(2, at = (1:nrow(l))/nrow(l), labels=rownames(l),tick=FALSE,las=2,cex.axis=0.2)
data(ants) allindivs <- c(union(as.character(ants$VertexFrom), as.character(ants$VertexTo)), "NULL1", "NULL2") l <- generatelatencies(ants, allindivs) image(l[,,1000],axes=FALSE,frame=TRUE,col=rainbow(100)) axis(1, at = (1:ncol(l))/ncol(l), labels=colnames(l),tick=FALSE,las=2,cex.axis=0.2) axis(2, at = (1:nrow(l))/nrow(l), labels=rownames(l),tick=FALSE,las=2,cex.axis=0.2)
Constructs weighted directed networks from all events occurring within certain time windows. Weight is equal to the number of interactions observed during the time window.
generatenetworkslices(g, timedeltas)
generatenetworkslices(g, timedeltas)
g |
The time-ordered network to be sliced. |
timedeltas |
A n x 2 matrix, where each row contains a set of start (first column) and stop (second column) times at which the network should be sliced. |
A list containing n time-aggregated networks corresponding to the n time windows.
Benjamin Blonder [email protected].
plotnetworkslices
, generatetimedeltas
, generatetimelags
~~~
data(ants) allindivs <- c(union(as.character(ants$VertexFrom), as.character(ants$VertexTo)), "NULL1", "NULL2") g <- generatetonetwork(ants, allindivs) td100 <- generatetimedeltas(0,1500,100) ns100 <- generatenetworkslices(g, td100) plotnetworkslices(ns100, td100)
data(ants) allindivs <- c(union(as.character(ants$VertexFrom), as.character(ants$VertexTo)), "NULL1", "NULL2") g <- generatetonetwork(ants, allindivs) td100 <- generatetimedeltas(0,1500,100) ns100 <- generatenetworkslices(g, td100) plotnetworkslices(ns100, td100)
-
generatetimeaggregatednetwork(g, starttime, stoptime)
generatetimeaggregatednetwork(g, starttime, stoptime)
g |
The time-ordered network to be aggregated |
starttime |
The time at which to begin aggregating interactions. |
stoptime |
The time at which to stop aggregating interactions. |
A weighted time-aggregated network whose edge weights equal the number of interactions between those vertices in the time window.
Benjamin Blonder [email protected].
data(ants) allindivs <- c(union(as.character(ants$VertexFrom), as.character(ants$VertexTo)), "NULL1", "NULL2") g <- generatetonetwork(ants, allindivs) tan500 <- generatetimeaggregatednetwork(g, 0, 500) plottanet(tan500)
data(ants) allindivs <- c(union(as.character(ants$VertexFrom), as.character(ants$VertexTo)), "NULL1", "NULL2") g <- generatetonetwork(ants, allindivs) tan500 <- generatetimeaggregatednetwork(g, 0, 500) plottanet(tan500)
-
generatetimedeltas(starttime, stoptime, delta)
generatetimedeltas(starttime, stoptime, delta)
starttime |
The starting time of the first time window. |
stoptime |
The stopping time of the last time window. |
delta |
The size of each time window. |
A n x 2 matrix. Each row contains the start and stop time of a window with width delta.
Benjamin Blonder [email protected].
generatetimelags
~~~
td100 <- generatetimedeltas(0,1500,100) boxplot(t(td100))
td100 <- generatetimedeltas(0,1500,100) boxplot(t(td100))
-
generatetimelags(starttime, stoptime, delta)
generatetimelags(starttime, stoptime, delta)
starttime |
The starting time of the first time window. |
stoptime |
The stopping time of the last time window. |
delta |
The size by which to increase each time window. |
A n x 2 matrix. Each row contains the start and stop time of a window with widths increasing by delta.
Benjamin Blonder [email protected].
tl100 <- generatetimelags(0,1500,100) boxplot(t(tl100))
tl100 <- generatetimelags(0,1500,100) boxplot(t(tl100))
Constructs a directed network describing the causally permitted paths between a set of vertices that interact at known times.
generatetonetwork(raw, allindivs)
generatetonetwork(raw, allindivs)
raw |
An event list, consisting of a data frame with four columns: VertexFrom, VertexTo, TimeStart, and TimeStop. Each row in this data frame represents a single directed interaction event between VertexFrom and VertexTo beginning at TimeStart and ending at TimeStop. |
allindivs |
A list of all possible vertices potentially including ones not observed interacting during the range of time reported in |
A weighted directed network of class 'igraph'. Each vertex represents an individual at a time during which an interaction occurred. Edges represent causally permitted paths of resource flow and have a TimeCost, describing the time between interactions for an individual, or is 0 if the edge represents an interaction, and a HopCost, which is 0 if the edge connects the same individual at multiple times and 1 if it connects different individuals at the same time.
Benjamin Blonder [email protected].
Kostakos V. Temporal Graphs. arXiv (2008) vol. physics.soc-ph
data(ants) allindivs <- c(union(as.character(ants$VertexFrom), as.character(ants$VertexTo)), "NULL1", "NULL2") g <- generatetonetwork(ants, allindivs) plottonet(g)
data(ants) allindivs <- c(union(as.character(ants$VertexFrom), as.character(ants$VertexTo)), "NULL1", "NULL2") g <- generatetonetwork(ants, allindivs) plottonet(g)
-
generatetonetworkfromvel(vel)
generatetonetworkfromvel(vel)
vel |
A data frame listing all directed edges |
Benjamin Blonder [email protected].
##---- Should be DIRECTLY executable !! ---- ##-- ==> Define data, use random, ##-- or do help(data=index) for the standard data sets.
##---- Should be DIRECTLY executable !! ---- ##-- ==> Define data, use random, ##-- or do help(data=index) for the standard data sets.
-
generatevertexedgelist(raw, allindivs)
generatevertexedgelist(raw, allindivs)
raw |
A data frame of events |
allindivs |
A vector of names |
Benjamin Blonder [email protected].
-
maxpoints(td)
maxpoints(td)
td |
A n x 2 matrix describing a set of start and stop times. |
A maximum value for each of n rows of td
Benjamin Blonder [email protected].
generatetimelags
,generatetimedeltas
~~~
tl100 <- generatetimelags(0,1500,100) boxplot(t(maxpoints(tl100)))
tl100 <- generatetimelags(0,1500,100) boxplot(t(maxpoints(tl100)))
-
midpoints(td)
midpoints(td)
td |
A n x 2 matrix describing a set of start and stop times. |
A mean value for each of n rows of td
Benjamin Blonder [email protected].
generatetimelags
,generatetimedeltas
~~~
tl100 <- generatetimelags(0,1500,100) boxplot(t(midpoints(tl100)))
tl100 <- generatetimelags(0,1500,100) boxplot(t(midpoints(tl100)))
-
plotnetworkslices(slices, timedeltas, ...)
plotnetworkslices(slices, timedeltas, ...)
slices |
A list of n time-aggregated networks |
timedeltas |
A n x 2 matrix describing the start and stop times for each time-aggregated network |
... |
Other arguments to be passed to |
None; used for its side effect of producing a plot.
Benjamin Blonder [email protected].
plotnetworkslices
, generatetimedeltas
, generatetimelags
~~~
data(ants) allindivs <- c(union(as.character(ants$VertexFrom), as.character(ants$VertexTo)), "NULL1", "NULL2") g <- generatetonetwork(ants, allindivs) td100 <- generatetimedeltas(0,1500,100) ns100 <- generatenetworkslices(g, td100) plotnetworkslices(ns100, td100)
data(ants) allindivs <- c(union(as.character(ants$VertexFrom), as.character(ants$VertexTo)), "NULL1", "NULL2") g <- generatetonetwork(ants, allindivs) td100 <- generatetimedeltas(0,1500,100) ns100 <- generatenetworkslices(g, td100) plotnetworkslices(ns100, td100)
Plots a time-aggregated network. See igraph.plotting for more details.
plottanet(timeaggregatednetwork, layout = layout.circle, vertex.label = V(timeaggregatednetwork)$name, vertex.size = 0, vertex.label.cex = 0.5, edge.arrow.size = 0.5, edge.width = E(timeaggregatednetwork)$Count/5, ...)
plottanet(timeaggregatednetwork, layout = layout.circle, vertex.label = V(timeaggregatednetwork)$name, vertex.size = 0, vertex.label.cex = 0.5, edge.arrow.size = 0.5, edge.width = E(timeaggregatednetwork)$Count/5, ...)
timeaggregatednetwork |
The network to print, an object of the igraph class |
layout |
Graph layout function - see ?layout in igraph |
vertex.label |
Vertex labels. Defaults to the name of each vertex. |
vertex.size |
Size of each vertex. |
vertex.label.cex |
Label size factor. |
edge.arrow.size |
Arrow size. |
edge.width |
Arrow width, defaults to be proportional to edge weight. |
... |
Other arguments to be passed to |
None; used for its side effect of producing a plot.
Benjamin Blonder [email protected].
data(ants) allindivs <- c(union(as.character(ants$VertexFrom), as.character(ants$VertexTo)), "NULL1", "NULL2") g <- generatetonetwork(ants, allindivs) tan <- generatetimeaggregatednetwork(g, 0, 500) plottanet(tan,layout=layout.kamada.kawai)
data(ants) allindivs <- c(union(as.character(ants$VertexFrom), as.character(ants$VertexTo)), "NULL1", "NULL2") g <- generatetonetwork(ants, allindivs) tan <- generatetimeaggregatednetwork(g, 0, 500) plottanet(tan,layout=layout.kamada.kawai)
Plots a time-ordered network with vertices ordinated along the x-axis and time increasing along the y-axis. Interactions are drawn as horizontal lines; vertices are connected to themselves in time by vertical lines.
plottonet(g, path = NULL, edgecolor = "gray", edgehighlightcolor = "red", vertex.size = 0.01, edge.arrow.size = 0.1, edge.width = 0.2, vertex.color = NA, vertex.label.cex = 0.1, vertex.frame.color = NA, vertex.label.color = "black")
plottonet(g, path = NULL, edgecolor = "gray", edgehighlightcolor = "red", vertex.size = 0.01, edge.arrow.size = 0.1, edge.width = 0.2, vertex.color = NA, vertex.label.cex = 0.1, vertex.frame.color = NA, vertex.label.color = "black")
g |
The time-ordered network to plot |
path |
If supplied, a particular list of vertices comprising a causally-permitted path that will be highlighted in the final illustration. |
edgecolor |
The color of all edges in the graph. |
edgehighlightcolor |
The color of the vertx path to be highlighted. |
vertex.size |
Vertex size. See igraph.plotting for more details. |
edge.arrow.size |
Edge arrow size. |
edge.width |
Edge width. |
vertex.color |
Vertex color. |
vertex.label.cex |
Vertex label size factor. |
vertex.frame.color |
Vertex frame color. |
vertex.label.color |
Vertex label color. |
None; used for its side-effect of producing a plot.
Benjamin Blonder [email protected].
data(ants) allindivs <- c(union(as.character(ants$VertexFrom), as.character(ants$VertexTo)), "NULL1", "NULL2") g <- generatetonetwork(ants, allindivs) plottonet(g)
data(ants) allindivs <- c(union(as.character(ants$VertexFrom), as.character(ants$VertexTo)), "NULL1", "NULL2") g <- generatetonetwork(ants, allindivs) plottonet(g)
edge_randomization
and randomized_edges
. An internal function.
NA
randomize_edges_helper(edges, randomize_vertices)
randomize_edges_helper(edges, randomize_vertices)
edges |
A data frame for an edge list |
randomize_vertices |
A binary variable |
Tim Gernat <[email protected]>
##---- Should be DIRECTLY executable !! ---- ##-- ==> Define data, use random, ##-- or do help(data=index) for the standard data sets. ## The function is currently defined as function (edges, randomize_vertices) { vertex_columns <- c("VertexFrom", "VertexTo") unique_edges <- unique(edges[, vertex_columns]) unique_edge_count <- nrow(unique_edges) edge_map <- cbind(unique_edges, unique_edges[sample(unique_edge_count, unique_edge_count), ]) new_vertex_columns <- c("NewVF", "NewVT") colnames(edge_map) <- c(vertex_columns, new_vertex_columns) if (randomize_vertices) { edge_map[, new_vertex_columns] <- sample(unlist(edge_map[, new_vertex_columns]), unique_edge_count * 2) repeat { invalid <- (edge_map$NewVF == edge_map$NewVT) | (duplicated(edge_map[, new_vertex_columns])) if (sum(invalid) == 0) break for (i in which(invalid)) edge_map <- swap(edge_map, i, sample(new_vertex_columns, 1), sample(unique_edge_count, 1), sample(new_vertex_columns, 1)) } } original_colnames <- colnames(edges) attribute_columns <- original_colnames[!(original_colnames %in% vertex_columns)] edges <- merge(edges, edge_map) edges <- edges[, c(new_vertex_columns, attribute_columns)] colnames(edges)[1:length(new_vertex_columns)] <- vertex_columns return(edges) }
##---- Should be DIRECTLY executable !! ---- ##-- ==> Define data, use random, ##-- or do help(data=index) for the standard data sets. ## The function is currently defined as function (edges, randomize_vertices) { vertex_columns <- c("VertexFrom", "VertexTo") unique_edges <- unique(edges[, vertex_columns]) unique_edge_count <- nrow(unique_edges) edge_map <- cbind(unique_edges, unique_edges[sample(unique_edge_count, unique_edge_count), ]) new_vertex_columns <- c("NewVF", "NewVT") colnames(edge_map) <- c(vertex_columns, new_vertex_columns) if (randomize_vertices) { edge_map[, new_vertex_columns] <- sample(unlist(edge_map[, new_vertex_columns]), unique_edge_count * 2) repeat { invalid <- (edge_map$NewVF == edge_map$NewVT) | (duplicated(edge_map[, new_vertex_columns])) if (sum(invalid) == 0) break for (i in which(invalid)) edge_map <- swap(edge_map, i, sample(new_vertex_columns, 1), sample(unique_edge_count, 1), sample(new_vertex_columns, 1)) } } original_colnames <- colnames(edges) attribute_columns <- original_colnames[!(original_colnames %in% vertex_columns)] edges <- merge(edges, edge_map) edges <- edges[, c(new_vertex_columns, attribute_columns)] colnames(edges)[1:length(new_vertex_columns)] <- vertex_columns return(edges) }
Produces a new event list from an existing event list with resampled vertex identities given certain constraints on randomization. Effectively re-orders pairs of From/To vertices between different times.
randomizeidentities(raw, withinvertexfrom, byvertexfrom, withreplacement)
randomizeidentities(raw, withinvertexfrom, byvertexfrom, withreplacement)
raw |
A raw event list to be resampled. Contains four columns: VertexFrom, VertexTo, TimeStart, TimeStop |
withinvertexfrom |
If true, resamples within data subsets where VertexFrom is fixed; otherwise resamples within all data. |
byvertexfrom |
If true, subsets of data for withinvertexfrom are obtained using VertexFrom; if false, using VertexTo. |
withreplacement |
Samples with or without replacement. |
An event list of the same size or smaller as raw. The returned event list will be smaller only if resampling produces events that connect a vertex to itself; these are removed.
Benjamin Blonder [email protected].
data(ants) allindivs <- c(union(as.character(ants$VertexFrom), as.character(ants$VertexTo)), "NULL1", "NULL2") ri <- randomizeidentities(ants,withinvertexfrom=TRUE,byvertexfrom=TRUE,withreplacement=TRUE) g <- generatetonetwork(ri, allindivs) plottonet(g)
data(ants) allindivs <- c(union(as.character(ants$VertexFrom), as.character(ants$VertexTo)), "NULL1", "NULL2") ri <- randomizeidentities(ants,withinvertexfrom=TRUE,byvertexfrom=TRUE,withreplacement=TRUE) g <- generatetonetwork(ri, allindivs) plottonet(g)
Produces a new event list from an existing event list with resampled event times given certain constraints on randomization. Effectively re-orders pairs of start/stop times between different vertices.
randomizetimes(raw, withinvertexfrom, byvertexfrom, withreplacement)
randomizetimes(raw, withinvertexfrom, byvertexfrom, withreplacement)
raw |
A raw event list to be resampled. Contains four columns: VertexFrom, VertexTo, TimeStart, TimeStop |
withinvertexfrom |
If true, resamples within data subsets where VertexFrom is fixed; otherwise resamples within all data. |
byvertexfrom |
If true, subsets of data for withinvertexfrom are obtained using VertexFrom; if false, using VertexTo. |
withreplacement |
Samples with or without replacement. |
An event list of the same size as raw with event times resampled. Resampling does not break the relationship between start and stop time; i.e. resampled events will have the same duration as original events.
Benjamin Blonder [email protected].
data(ants) allindivs <- c(union(as.character(ants$VertexFrom), as.character(ants$VertexTo)), "NULL1", "NULL2") rt <- randomizetimes(ants,withinvertexfrom=TRUE,byvertexfrom=TRUE,withreplacement=TRUE) g <- generatetonetwork(rt, allindivs) plottonet(g)
data(ants) allindivs <- c(union(as.character(ants$VertexFrom), as.character(ants$VertexTo)), "NULL1", "NULL2") rt <- randomizetimes(ants,withinvertexfrom=TRUE,byvertexfrom=TRUE,withreplacement=TRUE) g <- generatetonetwork(rt, allindivs) plottonet(g)
Take a data frame specifying the edges of a temporal network and create a randomized reference network which maintains certain properties of the original network and destroys others.
total_randomization(edges) randomly_permuted_times(edges) vertex_randomization(edges) contact_randomization(edges) time_reversal(edges) randomly_permuted_times(edges) random_times(edges) randomized_contacts(edges) edge_randomization(edges) randomized_edges(edges)
total_randomization(edges) randomly_permuted_times(edges) vertex_randomization(edges) contact_randomization(edges) time_reversal(edges) randomly_permuted_times(edges) random_times(edges) randomized_contacts(edges) edge_randomization(edges) randomized_edges(edges)
edges |
A |
randomly_permuted_times
permutes the start time of contacts and
adjusts the end time to maintain contact duration.
vertex_randomization
assigns vertices randomly and with equal
probability to contacts.
contact_randomization
randomly permutes vertices between contacts.
time_reversal
reverses the temporal order of contacts while
maintaining the temporal distance of contacts.
randomly_permuted_times
randomly permutes the start time of contacts
while maintaining contact duration.
random_times
assigns to the start time of each contact a random time
between min(edges$TimeStart)
and max(edges$TimeStop)
,
maintaining the duration of each contact.
randomized_contacts
redistributes contacts randomly among edges.
edge_randomization
randomly exchanges whole contact sequences between
edges.
randomized_edges
randomly rewires edges. When an edge gets rewired,
the contact sequence associated with that edge follow the edge.
total_randomization
assigns vertices randomly to contacts, assuming
that all vertices are equally likely participate in a contact
Randomized reference networks returned by these functions contain no contacts with self.
A data.frame
with the same columns as the edges
, specifying the
contacts of the randomized reference network.
Tim Gernat <[email protected]>
Holme & Saramaki, Physics Reports 519 (2012), p. 116-118
# load a temporal network require(timeordered) data(ants) # randomly permute contact start timestamps while preserving contact duration r1 <- randomly_permuted_times(ants) # randomly permute vertices between contacts and assign a random start # timestamp to each contact while preserving contact duration r2 <- contact_randomization(ants) r2 <- random_times(r1)
# load a temporal network require(timeordered) data(ants) # randomly permute contact start timestamps while preserving contact duration r1 <- randomly_permuted_times(ants) # randomly permute vertices between contacts and assign a random start # timestamp to each contact while preserving contact duration r2 <- contact_randomization(ants) r2 <- random_times(r1)
Randomly removes a fixed fraction of the event list.
rarefy(raw, fraction)
rarefy(raw, fraction)
raw |
The event list to be rarefied. |
fraction |
A fraction (between 0 and 1) of the events to be randomly deleted. |
An event list with floor(nrow(raw) * fraction) events remaining.
Benjamin Blonder [email protected].
randomizeidentities
,randomizetimes
##---- Should be DIRECTLY executable !! ---- ##-- ==> Define data, use random, ##-- or do help(data=index) for the standard data sets.
##---- Should be DIRECTLY executable !! ---- ##-- ==> Define data, use random, ##-- or do help(data=index) for the standard data sets.
-
shortesthoppath(g, startvertexname, startvertextime, stopvertexname, stopvertextime)
shortesthoppath(g, startvertexname, startvertextime, stopvertexname, stopvertextime)
g |
The time-ordered network on which to find paths. |
startvertexname |
The name of the start vertex. |
startvertextime |
The time of the start vertex. Must be a time at which an interaction has occurred involving this vertex. |
stopvertexname |
The name of the stop vertex. |
stopvertextime |
The time of the stop vertex. Must be a time at which an interaction has occurred involving this vertex. |
A vertex list containing all the events on the shortest-hop path between the start and stop vertices/times.
Multiple shortest-hop paths may exist; returns only one of them.
Benjamin Blonder [email protected].
data(ants) allindivs <- c(union(as.character(ants$VertexFrom), as.character(ants$VertexTo)), "NULL1", "NULL2") g <- generatetonetwork(ants, allindivs) shp <- shortesthoppath(g, "WBGG", 927, "GYGG", 1423) plottonet(g, shp) title(paste(length(unique(shp$Name))," hops"))
data(ants) allindivs <- c(union(as.character(ants$VertexFrom), as.character(ants$VertexTo)), "NULL1", "NULL2") g <- generatetonetwork(ants, allindivs) shp <- shortesthoppath(g, "WBGG", 927, "GYGG", 1423) plottonet(g, shp) title(paste(length(unique(shp$Name))," hops"))
-
shortesttimepath(g, startvertexname, startvertextime, stopvertexname)
shortesttimepath(g, startvertexname, startvertextime, stopvertexname)
g |
The time-ordered network on which to find paths. |
startvertexname |
The name of the start vertex. |
startvertextime |
The time of the start vertex. Must be a time at which an interaction has occurred involving this vertex. |
stopvertexname |
The name of the stop vertex. |
A vertex list containing all the events on the shortest-time path between the start vertex at the start time and the stop vertex at a later time.
May generate warning messages - don't worry!
Benjamin Blonder [email protected].
data(ants) allindivs <- c(union(as.character(ants$VertexFrom), as.character(ants$VertexTo)), "NULL1", "NULL2") g <- generatetonetwork(ants, allindivs) stp <- shortesttimepath(g, "WBGG", 927, "Q") plottonet(g, stp) title(paste(diff(range(stp$Time)), "time elapsed"))
data(ants) allindivs <- c(union(as.character(ants$VertexFrom), as.character(ants$VertexTo)), "NULL1", "NULL2") g <- generatetonetwork(ants, allindivs) stp <- shortesttimepath(g, "WBGG", 927, "Q") plottonet(g, stp) title(paste(diff(range(stp$Time)), "time elapsed"))
Determines the number of unique vertices that can be causally linked to an interaction event after a certain time delay. This function determines the fraction of unique vertices reached after a certain time from a random sample of interaction events.
spreadanalysis(g, timedelays, numsamples, normalizebyname=FALSE)
spreadanalysis(g, timedelays, numsamples, normalizebyname=FALSE)
g |
The time-ordered network to be studied. |
timedelays |
A vector time delays at which to determine the fraction of vertices reached. |
numsamples |
The number of random events to sample (without replacement) as seeds for the spreading process. |
normalizebyname |
If true, divides the number of vertices reached by the number of unique vertex names; if false, by the number of time-ordered vertices. |
A data frame whose columns are named for each time delay and contains the fraction of total vertices reached by a spreading process beginning from the seed vertices by the time delay.
Results can be aggregated by start vertex - see transformspreadbyindividual
Benjamin Blonder [email protected].
data(ants) allindivs <- c(union(as.character(ants$VertexFrom), as.character(ants$VertexTo)), "NULL1", "NULL2") g <- generatetonetwork(ants, allindivs) sa <- spreadanalysis(g, seq(0,1000,by=50), 20) boxplot(sa[,-1],xlab="Time delay",ylab="Fraction reached")
data(ants) allindivs <- c(union(as.character(ants$VertexFrom), as.character(ants$VertexTo)), "NULL1", "NULL2") g <- generatetonetwork(ants, allindivs) sa <- spreadanalysis(g, seq(0,1000,by=50), 20) boxplot(sa[,-1],xlab="Time delay",ylab="Fraction reached")
NA
swap(df, r1, c1, r2, c2)
swap(df, r1, c1, r2, c2)
df |
A dataframe |
r1 |
The first row to swap |
c1 |
The first column to swap |
r2 |
The second row to swap |
c2 |
The second column to swap |
Tim Gernat <[email protected]>
##---- Should be DIRECTLY executable !! ---- ##-- ==> Define data, use random, ##-- or do help(data=index) for the standard data sets. ## The function is currently defined as function (df, r1, c1, r2, c2) { tmp <- df[r1, c1] df[r1, c1] <- df[r2, c2] df[r2, c2] <- tmp return(df) }
##---- Should be DIRECTLY executable !! ---- ##-- ==> Define data, use random, ##-- or do help(data=index) for the standard data sets. ## The function is currently defined as function (df, r1, c1, r2, c2) { tmp <- df[r1, c1] df[r1, c1] <- df[r2, c2] df[r2, c2] <- tmp return(df) }
Converts a data frame of spreading samples into a data frame that is grouped by vertex identity.
transformspreadbyindividual(sa)
transformspreadbyindividual(sa)
sa |
A data frame returned by |
A data frame whose columns are the identities of vertices and whose rows are the mean fraction of vertices reached by the seed vertex at each time delay, averaged over all samples beginning at this vertex.
Benjamin Blonder [email protected].
data(ants) allindivs <- c(union(as.character(ants$VertexFrom), as.character(ants$VertexTo)), "NULL1", "NULL2") g <- generatetonetwork(ants, allindivs) sa <- spreadanalysis(g, seq(0,1000,by=50), 20) b <- transformspreadbyindividual(sa) plot(ts(b),plot.type="single",col=rainbow(ncol(b)),xlab="Time",ylab="Fraction reached") legend("bottomright",colnames(b),lwd=1,col=rainbow(ncol(b)),bg="white")
data(ants) allindivs <- c(union(as.character(ants$VertexFrom), as.character(ants$VertexTo)), "NULL1", "NULL2") g <- generatetonetwork(ants, allindivs) sa <- spreadanalysis(g, seq(0,1000,by=50), 20) b <- transformspreadbyindividual(sa) plot(ts(b),plot.type="single",col=rainbow(ncol(b)),xlab="Time",ylab="Fraction reached") legend("bottomright",colnames(b),lwd=1,col=rainbow(ncol(b)),bg="white")