#> Warning: package 'spam' was built under R version 4.0.5
#> Warning: package 'dotCall64' was built under R version 4.0.5
#> Warning: package 'RSQLite' was built under R version 4.0.5
#> Warning: package 'sp' was built under R version 4.0.5
Graphical abstract
set.seed(1)
OCN <- create_OCN(30, 20, outletPos = 1)
OCN <- aggregate_OCN(landscape_OCN(OCN), thrA = 3)
par(mfrow = c(1, 3), mai = c(0, 0, 0.2, 0.2))
draw_simple_OCN(OCN, thrADraw = 3)
title("Optimal Channel Network")
draw_elev3D_OCN(OCN, drawRiver = FALSE, addColorbar = FALSE, expand = 0.2, theta = -30)
title("Elevation")
draw_thematic_OCN(OCN$AG$streamOrder, OCN, discreteLevels = TRUE, colPalette = rainbow(4))
title("Strahler stream order")
OCNet enables the creation and analysis of Optimal Channel Networks (OCNs). These are oriented spanning trees (built on rectangular lattices made up of square pixels) that reproduce all scaling features characteristic of real, natural river networks (Rodriguez-Iturbe et al. 1992; Rinaldo et al. 2014). As such, they can be used in a variety of numerical and laboratory experiments in the fields of hydrology, ecology and epidemiology. Notable examples include studies on metapopulations and metacommunities (e.g. Carrara et al. 2012), scenarios of waterborne pathogen invasions (e.g. Gatto et al. 2013) and biogeochemichal processes in streams (e.g. Helton, Hall, and Bertuzzo 2018).
OCNs are obtained by minimization of a functional which represents total energy dissipated by water flowing through the network spanning the lattice. Such a formulation embeds the evidence that morphological and hydrological characteristics of rivers (in particular, water discharge and slope) follow a power-law scaling with drainage area. For an overview of the functionalities of the package, see Carraro et al. (2020). For details on the theoretical foundation of the OCN concept, see Rinaldo et al. (2014).
In graph theory, an oriented spanning tree is a subgraph of a graph \(G\) such that:
At the simplest aggregation level (flow direction - FD; see Section 4 below), OCNs are oriented spanning trees whose nodes are the pixels consituting the lattice and whose edges represent flow directions.
Moreover, OCNs, just like real rivers, are constituted of nodes whose indegree (i.e. the number of edges pointing towards a node) can assume any value while the outdegree (number of edges exiting from the node) is equal to 1, except for the root (or outlet node), whose outdegree is equal to 0. Nodes with null indegree are termed sources. Nodes with indegree larger than 1 are confluences.
OCNet also allows building multiple networks within a single lattice. Each of these networks is defined by its respective outlet, which represents the root of a subgraph; the union of all subgraphs contains all elements of \(G\). For simplicity, we will still refer to “OCNs” with regards to these multiple-outlet entities. In this case, strictly speaking, OCNs are not trees but rather forests.
An OCN is defined by an adjacency matrix \(\mathbf{W}\) with entries \(w_{ij}\) equal to 1 if node \(i\) drains into \(j\) and null otherwise. Owing to the previously described properties, all rows of \(\mathbf{W}\) have a single non-zero entry, except those identifying the outlet nodes, whose entries are all null. Each adjacency matrix uniquely defines a vector of contributing areas (or drainage areas) \(\mathbf{A}\), whose components \(A_i\) are equal to the number of nodes upstream of node \(i\) plus the node itself. Mathematically, this can be expressed as \((\mathbf{I}-\mathbf{W}^T)\mathbf{A}=\mathbf{1}\), where \(\mathbf{I}\) is the identity matrix and \(\mathbf{1}\) a vector of ones.
The generation of an OCN is performed by function create_OCN
. Its only required inputs are the dimensions of the rectangular lattice, but several other features can be implemented via its optional inputs (see the function documentation for details). The output of create_OCN
is a list, which can be used as input to the subsequent function landscape_OCN
, as shown by the dependency tree below (indentation of an item implies dependency on the function at the above level).
create_OCN
: performs the OCN search algorithm
landscape_OCN
: calculates the elevation field generated by the OCN
aggregate_OCN
: aggregates the OCN at various levels
draw_thematic_OCN
: draws OCN with colors of nodes depending on a themedraw_subcatchments_OCN
: draws partition of the lattice into subcatchments as a result of the aggregation process of the OCNpaths_OCN
: calculates paths and path lengths among OCN nodesrivergeometry_OCN
: evaluates hydraulic properties of the OCNOCN_to_igraph
: transforms the OCN into an igraph objectOCN_to_SSN
: transforms the OCN into a SpatialStreamNetwork objectdraw_contour_OCN
: draws “real-shaped”1 OCNs and catchment contoursdraw_elev2D_OCN
: draws 2D elevation fielddraw_elev3D_OCN
: draws 3D elevation fielddraw_elev3Drgl_OCN
: draws 3D elevation field (via rgl
rendering system)find_area_threshold_OCN
: finds relationship between threshold area and number of nodes at the RN and AG level (see Figure 4.1 for relevant definitions)draw_simple_OCN
: fast drawing of OCNscreate_peano
: creates Peano networksOCNet functions are intended to be applied in sequential order: for each non-drawing function, the input list is copied into the output list, to which new sublists and objects are added.
Adjacency matrices and contributing area vectors of an OCN can be defined at different aggregation levels. In the output of OCNet functions, variables characterizing the OCN at the different aggregation levels are grouped within separate sublists, each of them identified by a two-letter acronym (marked in bold in the list below). Figure 4.1 provides a graphical visualization of the correspondence among the main aggregation levels.
create_OCN
). Edges’ lengths are equal to either cellsize
(the size of a pixel side, optional input in create_OCN
) or cellsize*sqrt(2)
, depending on whether flow direction is horizontal/vertical or diagonal.A_thr
in aggregate_OCN
). Such a procedure is customary in the hydrological problem of extracting a river network based on digital elevation models of the terrain (O’Callaghan and Mark 1984), and corresponds to the geomorphological concept of erosion threshold (associated to a threshold in landscape-forming runoff, of which drainage area represents a proxy Rodriguez-Iturbe and Rinaldo 2001). Edges’ lengths are again equal to either cellsize
or cellsize*sqrt(2)
.Nodes at the AG level correspond to a subset of nodes at the RN level. In particular, nodes at the AG level belong to at least one of these four categories:
maxReachLength
is finite): nodes that split edges that are longer than maxReachLength
.Outlet nodes at the AG level might also be sources, confluences or breaking nodes. All AG nodes except outlet nodes have outdegree equal to 1. All RN nodes that do not correspond to AG nodes constitute the edges of the network at the AG level: more specifically, each edge is formed by an AG node and a sequence of RN nodes downstream of the AG node, until another AG node is found.
Figure 4.2 shows an alternative aggregation scheme for the network showed in Figure 4.1 when the optional input maxReachLength
is set to a finite value.
The output of aggregate_OCN
contains objects named OCN$XX$toYY
, where XX
and YY
are two different aggregation levels. These objects define the correspondences between indices among aggregation levels. OCN$XX$toYY
contains a number of elements equal to the number of nodes at XX
level; each element OCN$XX$toYY[[i]]
contains the index/indices at YY
level corresponding to node i
at XX
level. For aggregation level AG, additional correspondence objects are marked by the string Reach
: these consider the whole sequence of RN nodes constituting the edge departing from an AG node as belonging to the AG node.
The example shown in Figure 4.3 corresponds to the dataset OCN_4
included in the package. Note that index numbering starts from the lower-left (southwestern) corner of the lattice.
The R code below displays the different OCN$XX$toYY
objects corresponding to the example in Figure 4.3:
ex <- aggregate_OCN(landscape_OCN(OCN_4), thrA = 2)
ex$FD$toRN
#> [1] 1 2 3 0 4 0 0 5 6 7 0 0 0 0 8 0
ex$FD$toSC
#> [1] 1 3 3 3 2 3 3 3 4 5 5 3 4 5 5 5
ex$RN$toFD
#> [1] 1 2 3 5 8 9 10 15
ex$RN$toAG
#> [1] 1 0 0 2 3 4 0 5
ex$RN$toAGReach
#> [1] 1 3 3 2 3 4 5 5
ex$AG$toFD
#> [1] 1 5 8 9 15
ex$AG$ReachToFD
#> [[1]]
#> [1] 1
#>
#> [[2]]
#> [1] 5
#>
#> [[3]]
#> [1] 2 3 8
#>
#> [[4]]
#> [1] 9
#>
#> [[5]]
#> [1] 10 15
ex$AG$toRN
#> [1] 1 4 5 6 8
ex$AG$ReachToRN
#> [[1]]
#> [1] 1
#>
#> [[2]]
#> [1] 4
#>
#> [[3]]
#> [1] 2 3 5
#>
#> [[4]]
#> [1] 6
#>
#> [[5]]
#> [1] 7 8
ex$SC$toFD
#> [[1]]
#> [1] 1
#>
#> [[2]]
#> [1] 5
#>
#> [[3]]
#> [1] 2 3 8 4 6 7 12
#>
#> [[4]]
#> [1] 9 13
#>
#> [[5]]
#> [1] 10 15 11 14 16
Let’s build an OCN on a 20x20 lattice and assume that each cell represents a square of side 500 m. The total size of the catchment is therefore 100 km2. Let’s locate the outlet close to the southwestern corner of the lattice. Function draw_simple_OCN
can then be used to display the OCN.
set.seed(1)
OCNwe <- create_OCN(20, 20, outletPos = 3, cellsize = 500)
par(mai=c(0,0,0,0))
draw_simple_OCN(OCNwe)
Now, let’s construct the elevation field subsumed by the OCN. Let’s suppose that the outlet has null elevation and slope equal to 0.01. Then, we use draw_elev3D_OCN
to draw the three-dimensional elevation map (values are in m).
OCNwe <- landscape_OCN(OCNwe, slope0 = 0.01)
par(mai=c(0,0,0,0.5))
draw_elev3D_OCN(OCNwe, drawRiver = FALSE)
Next, the OCN can be aggregated. Let’s suppose that the desired number of nodes at the AG level be as close as possible4 to 20. With function find_area_threshold_OCN
we can derive the corresponding value of drainage area threshold:
thr <- find_area_threshold_OCN(OCNwe)
# find index corresponding to thr$Nnodes ~= 20
indThr <- which(abs(thr$nNodesAG - 20) == min(abs(thr$nNodesAG - 20)))
indThr <- max(indThr) # pick the last ind_thr that satisfies the condition above
thrA20 <- thr$thrValues[indThr] # corresponding threshold area
The resulting number of nodes is5 20, corresponding to a threshold area thrA20 =
2.5 km2. The latter value can now be used in function aggregate_OCN
to obtain the aggregated network. Function draw_subcatchments_OCN
shows how the lattice is partitioned into subcatchments. It is possible to add points at the locations of the nodes at the AG level.
OCNwe <- aggregate_OCN(OCNwe, thrA = thrA20)
par(mai=c(0.1,0,0.1,0))
draw_subcatchments_OCN(OCNwe)
points(OCNwe$AG$X,OCNwe$AG$Y, pch = 21, col = "blue", bg = "blue")
Finally, draw_thematic_OCN
can be used to display the along-stream distances of RN-level nodes to the outlet (in m), as calculated by paths_OCN
.
OCNwe <- paths_OCN(OCNwe, includePaths = TRUE)
#> RN downstream paths... 1.6%
RN downstream paths... 3.2%
RN downstream paths... 4.8%
RN downstream paths... 6.4%
RN downstream paths... 8.1%
RN downstream paths... 9.7%
RN downstream paths... 11.3%
RN downstream paths... 12.9%
RN downstream paths... 14.5%
RN downstream paths... 16.1%
RN downstream paths... 17.7%
RN downstream paths... 19.3%
RN downstream paths... 20.9%
RN downstream paths... 22.6%
RN downstream paths... 24.2%
RN downstream paths... 25.8%
RN downstream paths... 27.4%
RN downstream paths... 29.0%
RN downstream paths... 30.6%
RN downstream paths... 32.2%
RN downstream paths... 33.8%
RN downstream paths... 35.4%
RN downstream paths... 37.1%
RN downstream paths... 38.7%
RN downstream paths... 40.3%
RN downstream paths... 41.9%
RN downstream paths... 43.5%
RN downstream paths... 45.1%
RN downstream paths... 46.7%
RN downstream paths... 48.3%
RN downstream paths... 50.0%
RN downstream paths... 51.6%
RN downstream paths... 53.2%
RN downstream paths... 54.8%
RN downstream paths... 56.4%
RN downstream paths... 58.0%
RN downstream paths... 59.6%
RN downstream paths... 61.2%
RN downstream paths... 62.8%
RN downstream paths... 64.5%
RN downstream paths... 66.1%
RN downstream paths... 67.7%
RN downstream paths... 69.3%
RN downstream paths... 70.9%
RN downstream paths... 72.5%
RN downstream paths... 74.1%
RN downstream paths... 75.7%
RN downstream paths... 77.3%
RN downstream paths... 79.0%
RN downstream paths... 80.6%
RN downstream paths... 82.2%
RN downstream paths... 83.8%
RN downstream paths... 85.4%
RN downstream paths... 87.0%
RN downstream paths... 88.6%
RN downstream paths... 90.2%
RN downstream paths... 91.8%
RN downstream paths... 93.5%
RN downstream paths... 95.1%
RN downstream paths... 96.7%
RN downstream paths... 98.3%
RN downstream paths... 99.9%
RN downstream paths... 100.0%
#> AG downstream paths... 5.0%
AG downstream paths... 10.0%
AG downstream paths... 15.0%
AG downstream paths... 20.0%
AG downstream paths... 25.0%
AG downstream paths... 30.0%
AG downstream paths... 35.0%
AG downstream paths... 40.0%
AG downstream paths... 45.0%
AG downstream paths... 50.0%
AG downstream paths... 54.9%
AG downstream paths... 59.9%
AG downstream paths... 64.9%
AG downstream paths... 69.9%
AG downstream paths... 74.9%
AG downstream paths... 79.9%
AG downstream paths... 84.9%
AG downstream paths... 89.9%
AG downstream paths... 94.9%
AG downstream paths... 99.9%
AG downstream paths... 100.0%
par(mai=c(0.1,0,0.1,0))
draw_thematic_OCN(OCNwe$RN$downstreamPathLength[ , OCNwe$RN$outlet], OCNwe,
backgroundColor = "#606060")
Let’s build a simple discrete-time, deterministic metapopulation model on the previously built OCN. In particular, let’s assume that:
Therefore, the model equation is: \[ \begin{split} P_i(t+1) &= \frac{r P_i(t)}{1+ (r-1)P_i(t)/K_i} - (p_d D_i + p_u U_i)G(K_i,K_i)\frac{P_i(t)}{K_i} \\ &\quad + p_d \left( \sum_{j=1} w_{ji} G(K_j,K_j)\frac{P_j(t)}{K_j}\right) + p_u \left( \sum_{j=1} w_{ij} Y_i G(K_j,K_j)\frac{P_j(t)}{K_j}\right) \end{split} \] where \(D_i\) (\(U_i\)) is equal to one if there is a downstream (upstream) connection available from node \(i\) and is null otherwise. Weights \(Y_i\) are defined as: \[ Y_i = \frac{A_i}{\sum_{k=1} w_{kj}A_k}, \] where \(j\) identifies the node downstream of \(i\). Moreover, it is \(Y_o = 1\).
At carrying capacity, the system is at equilibrium, which implies that the (expected) number of individuals moving from a node \(i\) to its downstream connection \(j\) is equal to the (expected) number of individuals moving from \(j\) to \(i\): \[ p_d G(K_i,K_i) = p_u Y_i G(K_j,K_j). \] The iterative application of the above equation allows the calculation of \(G(K_i,K_i)\) for all \(i\) up to a constant. To this end, let’s assume \(G(K_o,K_o) = G_o\). We therefore obtain \[ G(K_i,K_i) = G_o \left(\frac{p_u}{p_d}\right)^{|P_{io}|} \prod_{k \in P_{io}} Y_k, \] where \(P_{io}\) is the set of nodes constituting the downstream path from \(i\) to the outlet \(o\), while \(|P_{io}|\) is its cardinality.
For this example, let’s use the previously derived OCN_we
aggregated at the RN level. Let’s assume that carrying capacity is proportional to the river width evaluated at the nodes, while the initial population distribution is randomly assigned.
## Input data
OCNwe <- rivergeometry_OCN(OCNwe, widthMax = 5) # evaluate river width
K <- 10*OCNwe$RN$width # calculate carrying capacity
pop0 <- 2*mean(K)*runif(OCNwe$RN$nNodes) # initial random population vector
nTimestep <- 100 # number of timesteps
r <- 1.05 # proliferation rate
pd <- 0.5 # probability to move downstream
pu <- 1 - pd # probability to move upstream
Go <- 5 # parameter controlling mobility
# (no. individuals exiting from outlet node at carrying capacity is pu*Go)
We can now compute weights \(Y\):
## Weights for upstream movement
Y <- rep(1,OCNwe$RN$nNodes)
for (i in 1:OCNwe$RN$nNodes){
if (i != OCNwe$RN$outlet){
Y[i] <- OCNwe$RN$A[i]/(OCNwe$RN$W[ , OCNwe$RN$downNode[i]] %*% OCNwe$RN$A)
}
}
and \(G(K_i,K_i)\):
## Evaluate expected number of individuals moving at carrying capacity
GKK <- rep(0, OCNwe$RN$nNodes)
for (i in (1:OCNwe$RN$nNodes)){
path <- OCNwe$RN$downstreamPath[[i]][[OCNwe$RN$outlet]] # select path
prod <- Go # initialize product of Y
for (j in path){
prod <- prod*Y[j]
}
GKK[i] <- (pu/pd)^(length(path))*prod
}
We can now run the metapopulation model:
## Run metapopulation model
pop <- matrix(data=0,ncol=nTimestep,nrow=OCNwe$RN$nNodes) # metapopulation matrix
pop[,1] <- pop0 # initialization
for (t in 2:nTimestep){
for (i in 1:OCNwe$RN$nNodes){
pop[i, t] <-
# Beverton-Holt growth model
r*pop[i, t-1]/(1 + pop[i, t-1]*(r-1)/K[i]) +
# individuals exiting from node i
- (pu*(sum(OCNwe$RN$W[ , i])>0) + pd*(sum(OCNwe$RN$W[i, ])>0)) *
GKK[i] * (pop[i,t-1]/K[i]) +
# individuals entering in i from the upstream nodes
+ pd * OCNwe$RN$W[ , i] %*% (GKK*pop[ , t-1]/K) +
# individuals entering in i from the downstream node
+ pu * Y[i] * OCNwe$RN$W[i, ] %*% (GKK*pop[ , t-1]/K)
}
}
The left panel of the graph below shows the time evolution of the local population at the outlet (red) and at the node at highest distance from the outlet (blue). In the right panel, the evolution total metapopulation size is shown. Dashed lines indicate population values at carrying capacity.
par(mfrow = c(1, 2))
plot(pop[OCNwe$RN$outlet, ], type = "l", ylim = c(0, 1.05*K[OCNwe$RN$outlet]), col = "red",
xlab = "Time", ylab = "Population", lwd = 2)
title("Evolution of local pop. size")
lines(c(1, nTimestep),c(K[OCNwe$RN$outlet], K[OCNwe$RN$outlet]), col = "red", lty = 2)
farthestNode <- which(OCNwe$RN$downstreamPathLength[ , OCNwe$RN$outlet]
== max(OCNwe$RN$downstreamPathLength[ , OCNwe$RN$outlet]))[1]
lines(pop[farthestNode, ], type="l", col="blue",lwd=2)
lines(c(1, nTimestep), c(K[farthestNode], K[farthestNode]), col = "blue", lty = 2)
plot(colSums(pop), type = "l", xlab = "Time", ylab = "Population", lwd = 2, ylim = c(0, 1.05*sum(K)))
lines(c(1, nTimestep), c(sum(K),sum(K)), lty = 2)
title("Evolution of metapop. size")
Function draw_thematic_OCN
can be used to visualize the spatial distribution of the metapopulation at given time points.
par(mfrow = c(2, 2), mai = c(0.1, 0, 0.2, 0))
draw_thematic_OCN(pop[,1], OCNwe, colLevels = c(0, max(K), 1000),
drawNodes = TRUE)
title("Time = 1")
draw_thematic_OCN(pop[,5], OCNwe, colLevels = c(0, max(K), 1000),
drawNodes = TRUE)
title("Time = 5")
draw_thematic_OCN(pop[,20], OCNwe, colLevels = c(0, max(K), 1000),
drawNodes = TRUE)
title("Time = 20")
draw_thematic_OCN(pop[,100], OCNwe, colLevels = c(0, max(K), 1000),
drawNodes = TRUE)
title("Time = 100")
Function create_peano
can be used in lieu of create_OCN
to generate Peano networks on square lattices. Peano networks are deterministic, plane-filling fractals whose topological measures (Horton’s bifurcation and length ratios) are akin to those of real river networks (Marani, Rigon, and Rinaldo 1991) and can then be used in a variety of synthetic experiments, as it is the case for OCNs (e.g. Campos, Fort, and Méndez 2006). Peano networks are generated by means of an iterative algorithm: at each iteration, the size of the lattice side is doubled (see code below). As a result, Peano networks span squares of side equal to a power of 2. The outlet must be located at a corner of the square.
par(mfrow = c(2, 3), mai = c(0, 0, 0.2, 0))
peano0 <- create_peano(0)
draw_simple_OCN(peano0)
title("Iteration: 0 - Lattice size: 2x2")
peano1 <- create_peano(1)
draw_simple_OCN(peano1)
title("Iteration: 1 - Lattice size: 4x4")
peano2 <- create_peano(2)
draw_simple_OCN(peano2)
title("Iteration: 2 - Lattice size: 8x8")
peano3 <- create_peano(3)
draw_simple_OCN(peano3)
title("Iteration: 3 - Lattice size: 16x16")
peano4 <- create_peano(4)
draw_simple_OCN(peano4)
title("Iteration: 4 - Lattice size: 32x32")
peano5 <- create_peano(5)
draw_simple_OCN(peano5)
title("Iteration: 5 - Lattice size: 64x64")
The output of create_peano
is a list containing the same objects as those produced by create_OCN
. As such, it can be used as input for all other complementary functions of the package.
OCNet contains some ready-made large OCNs built via function create_OCN
. Their features are summarized in the Table below. Refer to the documentation of create_OCN
for the definition of column names. Note that:
coolingRate = 10
, initialNoCoolingPhase = 0
;coolingRate = 1
, initialNoCoolingPhase = 0
;coolingRate = 0.5
, initialNoCoolingPhase = 0.1
.set.seed
prior to executing create_OCN
.name | dimX | dimY | No. of outlets | Periodic Boundaries | Initial State | Cooling Schedule | Cellsize | seed | On CRAN? |
---|---|---|---|---|---|---|---|---|---|
OCN_4 7 |
4 | 4 | 1 | FALSE | 1 | ||||
OCN_20 |
20 | 20 | 1 | FALSE | I | default | 1 | 1 | Yes |
OCN_250 |
250 | 250 | 1 | FALSE | I | default | 1 | 2 | No |
OCN_250_T |
250 | 250 | 1 | FALSE | T | default | 1 | 2 | Yes |
OCN_250_V |
250 | 250 | 1 | FALSE | V | default | 1 | 2 | No |
OCN_250_cold |
250 | 250 | 1 | FALSE | I | cold | 1 | 2 | No |
OCN_250_hot |
250 | 250 | 1 | FALSE | I | hot | 1 | 2 | No |
OCN_250_V_cold |
250 | 250 | 1 | FALSE | V | cold | 1 | 2 | No |
OCN_250_V_hot |
250 | 250 | 1 | FALSE | V | hot | 1 | 2 | No |
OCN_250_PB |
250 | 250 | 1 | TRUE | I | default | 1 | 2 | Yes |
OCN_rect1 |
450 | 150 | 1 | FALSE | I | default | 1 | 3 | No |
OCN_rect2 |
150 | 450 | 1 | FALSE | I | default | 1 | 3 | No |
OCN_300_diag |
300 | 300 | 18 | FALSE | V | default | 50 | 4 | No |
OCN_300_4out |
300 | 300 | 49 | FALSE | V | default | 50 | 5 | Yes |
OCN_300_4out_PB_hot |
300 | 300 | 410 | TRUE | V | hot | 50 | 5 | Yes |
OCN_300_7out |
300 | 300 | 7 | FALSE | V | default | 50 | 5 | No |
OCN_400_T_out |
400 | 400 | 1 | FALSE | T | hot | 50 | 7 | No |
OCN_400_Allout |
400 | 400 | All | FALSE | H | hot | 50 | 8 | Yes |
OCN_500_hot 11 |
500 | 500 | 112 | FALSE | I | hot | 50 | 9 | No |
OCN_500_PB_hot |
500 | 500 | 1 | TRUE | V | hot | 50 | 10 | No |
Adjacency matrices at all aggregation levels are produced as spam
(Furrer and Sain 2010) objects. In order to transform the OCN into an igraph
(Csardi and Nepusz 2006) graph object, the adjacency matrix must be converted into a Matrix
object (via function as.dgCMatrix.spam
of spam
). Function graph_from_adjacency_matrix
of igraph
can then be used to obtain a graph object.
For example, let’s transform the previously obtained OCN_we
at the AG level into a graph:
par(mai=c(0.1,0.1,0.1,0.1))
g <- OCN_to_igraph(OCNwe, level = "AG")
plot.igraph(g, vertex.color = rainbow(OCNwe$AG$nNodes),
layout = matrix(c(OCNwe$AG$X,OCNwe$AG$Y),ncol = 2, nrow = OCNwe$AG$nNodes))
The same network can be displayed as an OCN:
Function OCN_to_SSN
transforms an OCN at a given aggregation level into an SSN (Ver Hoef et al. 2014) object. See the following example:
ssnOCN <- OCN_to_SSN(OCNwe, level = "RN", obsDesign = SSN::binomialDesign(50),
path = paste(tempdir(), "/OCN.ssn", sep = ""), importToR = TRUE)
plot.SpatialStreamNetwork(ssnOCN, "upDist", breaktype = "user", brks = seq(0,14000,2000),
xlab = "x [m]", ylab = "y [m]", asp = 1)
title("Distance from outlet of observation points [m]")
Campos, D., J. Fort, and V. Méndez. 2006. “Transport on Fractal River Networks. Application to Migration Fronts.” Theoretical Population Biology 69(1): 88–93. http://dx.doi.org/10.1016/j.tpb.2005.09.001.
Carrara, F., F. Altermatt, I. Rodriguez-Iturbe, and A. Rinaldo. 2012. “Dendritic Connectivity Controls Biodiversity Patterns in Experimental Metacommunities.” Proceedings of the National Academy of Sciences of the United States of America 109(15): 5761–6. http://dx.doi.org/10.1073/pnas.1119651109.
Carraro, L., E. Bertuzzo, E. A. Fronhofer, R. Furrer, I. Gounand, A. Rinaldo, and F. Altermatt. 2020. “Generation and Application of River Network Analogues for Use in Ecology and Evolution.” bioRxiv. http://dx.doi.org/10.1101/2020.02.17.948851.
Csardi, G., and T. Nepusz. 2006. “The Igraph Software Package for Complex Network Research.” InterJournal Complex Systems: 1695. https://igraph.org.
Furrer, R., and S. R. Sain. 2010. “Spam. A Sparse Matrix R Package with Emphasis on MCMC Methods for Gaussian Markov Random Fields.” Journal of Statistical Software 36(10): 1–25. http://dx.doi.org/10.18637/jss.v036.i10.
Gatto, M., L. Mari, E. Bertuzzo, R. Casagrandi, L. Righetto, I. Rodriguez-Iturbe, and A. Rinaldo. 2013. “Spatially Explicit Conditions for Waterborne Pathogen Invasion.” American Naturalist 182(3): 328–46. http://dx.doi.org/10.1086/671258.
Helton, A. M., R. O. Hall, and E. Bertuzzo. 2018. “How Network Structure Can Affect Nitrogen Removal by Streams.” Freshwater Biology 63(1): 128–40. http://dx.doi.org/10.1111/fwb.12990.
Marani, A., R. Rigon, and A. Rinaldo. 1991. “A Note on Fractal Channel Networks.” Water Resources Research 27(12): 3041–9. http://dx.doi.org/10.1029/91WR02077.
O’Callaghan, J. F., and D. A. Mark. 1984. “The Extraction of the Drainage Networks from Digital Elevation Data.” Computer Vision, Graphics, and Image Processing 28: 323–44. http://dx.doi.org/10.1016/S0734-189X(84)80011-0.
Rinaldo, A., R. Rigon, J. R. Banavar, A. Maritan, and I. Rodriguez-Iturbe. 2014. “Evolution and Selection of River Networks. Statics, Dynamics, and Complexity.” Proceedings of the National Academy of Sciences of the United States of America 111(7): 2417–24. http://dx.doi.org/10.1073/pnas.1322700111.
Rodriguez-Iturbe, I., and A. Rinaldo. 2001. Fractal River Basins. Chance and Self-Organization. Cambridge University Press.
Rodriguez-Iturbe, I., A. Rinaldo, R. Rigon, R. L. Bras, A. Marani, and E. Ijjász-Vásquez. 1992. “Energy Dissipation, Runoff Production, and the Three‐dimensional Structure of River Basins.” Water Resources Research 28(4): 1095–1103. http://dx.doi.org/10.1029/91WR03034.
Ver Hoef, J. M., E. E. Peterson, D. Clifford, and R. Shah. 2014. “SSN. An R Package for Spatial Statistical Modeling on Stream Networks.” Journal of Statistical Software 56(3): 1–45. http://dx.doi.org/10.18637/jss.v056.i03.
If the OCN was generated with option periodicBoundaries = FALSE
, the real-shaped representation of the OCN is equal to the simple, within-lattice representation provided by draw_simple_OCN
. If the OCN was generated with periodicBoundaries = TRUE
, its real-shaped representation is such that, starting from the outlet(s), the position of a FD node i
is flipped with respect to a lattice side if the edge between i
and the pixel downstream of i
crosses such lattice side.↩︎
This might not be the case when there are many outlets. If the area drained by one of these outlets is lower than the threshold imposed, the cluster of pixels drained by that outlet constitutes a node at the SC level that does not have a correspondence at the AG level.↩︎
Note that this pattern of flow directions was not derived by an OCN search algortihm, but rather drawn manually for illustration purposes.↩︎
It is not guaranteed that an OCN can be aggregated to an arbitrary number of nodes. This is due to the fact that the indegree of confluence nodes is typically 2 or 3, depending on the flow direction patterns. For example, if the outlet is not a confluence node, all confluence nodes have indegree equal to 2, and maxReachLength = Inf
, then the resulting number of aggregated nodes will be even.↩︎
The number of aggregated nodes is here uneven because one node has indegree equal to 3 (see figure below).↩︎
At this GitHub link it is possible to download the same OCNs generated by adding optional arguments saveEnergy = TRUE
, saveExitFlag = TRUE
, which substantially increases the size of the output of create_OCN
.↩︎
This is actually not an OCN, but was rather generated manually for illustration purposes.↩︎
outletPos = 1
.↩︎
outletSide = c("S","N","W","E")
, outletPos = c(1,300,149,150)
.↩︎
outletSide = c("S","N","W","E")
, outletPos = c(1,300,149,150)
.↩︎
This OCN was obtained with nIter = 50*dimX*dimY
.↩︎
outletSide = "E", outletPos = 100
↩︎