======================================================
======================================================
Since the spdep package was created, spatial weights objects have been constructed as lists with three components and a few attributes, in old-style class listw
objects. The first component of a listw
object is an nb
object, a list of n
integer vectors, with at least a character vector region.id
attribute with n
unique values (like the row.names
of a data.frame
object); n
is the number of spatial entities. Component i
of this list contains the integer identifiers of the neighbours of i
as a sorted vector with no duplication and values in 1:n
; if i
has no neighbours, the component is a vector of length 1
with value 0L
. The nb
object may contain an attribute indicating whether it is symmetric or not, that is whether i
is a neighbour of j
implies that j
is a neighbour of i
. Some neighbour definitions are symmetric by construction, such as contiguities or distance thresholds, others are asymmetric, such as k
-nearest neighbours. The nb
object redundantly stores both i
-j
and j
-i
links.
The second component of a listw
object is a list of n
numeric vectors, each of the same length as the corresponding non-zero vectors in the nb
object. These give the values of the spatial weights for each i
-j
neighbour pair. It is often the case that while the neighbours are symmetric by construction, the weights are not, as for example when weights are row-standardised by dividing each row of input weights by the count of neighbours or cardinality of the neighbour set of i
. In the nb2listw
function, it is also possible to pass through general weights, such as inverse distances, shares of boundary lengths and so on.
The third component of a listw
object records the style
of the weights as a character code, with "B"
for binary weights taking values zero or one (only one is recorded), "W"
for row-standardised weights, and so on. In order to subset listw
objects, knowledge of the style
may be necessary
It is obvious that this is similar to the way in which sparse matrices are stored, either by row - like the listw
object, or by column. The key insight is that storing zero values is unnecessary, as we only need to store the row and column locations of non-zero values. Early on, a Netlib library was used to provide limited support in spdep for sparse matrices, followed by functionality in SparseM, spam, and Matrix.
Since Matrix is a recommended package, its functionality has increasingly been used over time, and it has become one of two packages on which spdep depends. This is reported on loading:
library(spdep)
The legacy Columbus OH data set has 49 spatial entities, polygons, defined as the boundaries of policing districts in the city. spdep suggests maptools for portability, but in regular applications, rgdal should be used, because it handles coordinate reference systems more gracefully. We see that the shapefile-based entity IDs are zero-based:
library(maptools) columbus <- readShapePoly(system.file("etc/shapes/columbus.shp", package="spdep")[1]) row.names(columbus)[1:10]
Contiguous neighbours are often used for polygonal spatial entities, here with the poly2nb function defaulting to the queen criterion - entities are neighbours if they share a boundary point. We see that the entity IDs are copied across to the nb
object:
nb_q <- poly2nb(columbus) nb_q attr(nb_q, "region.id")[1:10] is.symmetric.nb(nb_q)
In order to make the object more complicated, let us drop the neighbour links for the 21st entity (noting that the print method reports the ID of the entity with no neighbours, not its number in 1:n
), and plot the resulting map of neighbours:
col2 <- droplinks(nb_q, 21) nb_q[[21]] col2[[21]] col2 is.symmetric.nb(col2) coords <- coordinates(columbus) plot(nb_q, coords, col="grey") plot(col2, coords, add=TRUE)
======================================================
At present only listw
objects can be coerced to objects of classes defined in Matrix. Because the style
is lost on coercion, it may not be possible to reconstruct spatial weights as the sparse matrix representation does not preserve it. We will start with symmetric binary weights, first creating a spatial weights object, and signalling that one entity has no neighbours with the zero.policy
argument (default false). The matrix and graph representations of no-neighbour entities are not obvious.
nb_B <- nb2listw(col2, style="B", zero.policy=TRUE) nb_B$style
spdep provides coercion methods from listw
to the "symmetricMatrix"
, "RsparseMatrix"
and "CsparseMatrix"
classes defined in Matrix. The "RsparseMatrix"
is the representation that is most similar to listw
, as it is row-based, but it is used less frequently in operations on sparse matrices. The entity IDs are passed using sparse matrix row and column names at present. Here we believe that our listw
object can be represented as a symmetric matrix, storing only a triangle rather than both i
-j
and j
-i
weights. The coercion method does check whether symmetry is present before proceeding:
B <- as(nb_B, "symmetricMatrix") all(B == t(B)) str(B) rownames(B)[1:10]
Let us now try to retreive the list of neighbours from the symmetric sparse matrix. At present, we have to coerce from one Matrix internal representation to another in order to get to the "dgCMatrix"
format used inside mat2listw
, so we coerce to "dgTMatrix"
from "dsTMatrix"
. The style is not retreived automatically, but is set to "M"
to indicate conversion from a matrix. The neighbour links are retreived correctly, as are the IDs:
nb_B1 <- mat2listw(as(B, "dgTMatrix")) nb_B1$style all.equal(nb_B1$neighbours, col2, check.attributes=FALSE) all.equal(attr(nb_B1$neighbours, "region.id"), attr(nb_B$neighbours, "region.id"))
An initial reason for implementing support for sparse weights matrices in spdep was to permit the calculation of the log determinant term in spatial regressions for larger data sets. Using the eigenvalue approach with for example eigenw
is limited by the need to operate on dense matrices in memory to solve the eigenproblem:
rho <- 0.1 sum(log(1 - rho * eigenw(nb_B)))
When n
is large, this may become impractical and/or time-consuming, but does permit the rapid calculation of values of the log determinant for differing values of the spatial coefficient ( \rho ). The Matrix package provides many determinant
methods, here for a "dsCMatrix"
resulting from subtracting a "dsCMatrix"
, the product of a scalar and a "dsTMatrix"
, from a "ddiMatrix"
. The value of the log determinant follows, calling a sparse Cholesky decomposition internally for suitable input matrices.
n <- nrow(B) I <- Diagonal(n) class(I - rho * B) c(determinant(I - rho * B, logarithm=TRUE)$modulus)
The computation of a sparse Cholesky decomposition for each value of the spatial coefficient ( \rho ) may be avoided by updating a pre-computed object; this approach provides fast and accurate log determinants for larger (but not very large) data sets:
nW <- -B nChol <- Cholesky(nW, Imult=8) n * log(rho) + (2 * c(determinant(update(nChol, nW, 1/rho))$modulus))
The use of row-standardisation leads to asymmetry even if the underlying neighbours are symmetric, unless all entities have matching numbers of neighbours (for example a regular grid on a torus):
nb_W <- nb2listw(col2, style="W", zero.policy=TRUE) W <- as(nb_W, "CsparseMatrix") str(W) all(W == t(W))
The lag
method for listw
objects is often used to create spatially lagged values, and returns the same values as the vector given by the product of the sparse general matrix and an input numeric vector. Note that by setting zero.policy
to TRUE
, the spatial lag of entity 21, which has no neighbours, is zero by construction:
set.seed(1) x <- runif(n) r1 <- as.numeric(W %*% x) r2 <- lag(nb_W, x, zero.policy=TRUE) all.equal(r1, r2, check.attributes=FALSE) plot(x, r1, ylim=c(0,1)) c(x[21], r1[21])
Calculating the log determinant for asymmetric weights (here with symmetric neighbours and symmetry induced by non-constant numbers of neighbours) may be carried out using eigenvalues as before, but the result may be a complex vector (here it is not, as discussed below). The appropriate determinant
method for "dgCMatrix"
objects uses an LU decomposition internally:
rho <- 0.5 sum(log(1 - rho * eigenw(nb_W))) class(I - rho * W) c(determinant(I - rho * W, logarithm=TRUE)$modulus)
We can show the internal workings of the method as:
LU <- lu(I - rho * W) sum(log(abs(diag(slot(LU, "U")))))
The nb2listw
function stores components that can be employed to transform the asymmetric weights matrix to symmetry by similarity, permitting the same log determinant to be computed using less costly numerical methods. The "W"
style used the cardinalities of neighbour sets (row sums) to introduce row standardisation, and they are stored as an attribute:
d <- attr(nb_W$weights, "comp")$d all.equal(d, card(col2))
If we first restore the row-standarised matrix to its binary form (which must be symmetric), we can pre- and post-multiply by the square roots of the inverted neighbour counts, yielding a symmetric matrix with the appropriate properties:
dW <- Diagonal(n, d) %*% W all(dW == t(dW)) isd <- Diagonal(n, 1/sqrt(d)) isd[21,21] Ws <- as(isd %*% dW %*% isd, "symmetricMatrix") rowSums(Ws)[21] class(Ws) c(determinant(I - rho * Ws, logarithm=TRUE)$modulus)
As can be seen, the division by the square root of zero for entity 21 is not a problem as the row of dW
is zero. The transformation by similarity permits the use of numerical methods for sparse symmetric matrices (and equivalently for eigenvalues and dense matrices). Note that this transformation is not available for intrinsically asymmetric neighbours, or for intrinsically asymmetric general weights.
arpack
in igraph for finding some eigenvaluesIn spatial regression, the domain of the spatial coefficient is given by the inverse of the maximum and minimum eigenvalues. When n
is moderate, we have the eigenvalues anyway, so the interval for line search is available without extra effort. When n
is somewhat larger, use may be made of the arpack
function in igraph. The function is not stable and requires manual tuning, so its use is commented out here:
1/range(eigenw(nb_B)) library(igraph) #f2 <- function(x, extra=NULL) {as.vector(B %*% x)} #1/arpack(f2, sym=TRUE, options=list(n=n, nev=1, ncv=8, which="LA", maxiter=200))$values #1/arpack(f2, sym=TRUE, options=list(n=n, nev=1, ncv=8, which="SA", maxiter=200))$values # At line 409 of file dsaupd.f: Fortran runtime error: Actual string # length is shorter than the declared one for dummy argument 'which' (0/2) #1/arpack(f2, sym=TRUE, options=list(n=n, nev=2, ncv=8, which="BE", maxiter=200))$values # "BE" gives: At line 558 of file dsaup2.f: Fortran runtime error: # Index '9' of dimension 1 of array 'bounds' above upper bound of 8
In this case, the results are trivial with small n
and a symmetric sparse matrix. It can be challenging to find the maxiter=
cutoff for larger n
; the which="BE"
argument takes both ends here.
1/range(eigenw(nb_W)) f2 <- function(x, extra=NULL) {as.vector(W %*% x)} #1/arpack(f2, sym=FALSE, options=list(n=n, nev=1, ncv=8, which="LR", maxiter=200))$values #1/arpack(f2, sym=FALSE, options=list(n=n, nev=1, ncv=8, which="SR", maxiter=200))$values # At line 409 of file dsaupd.f: Fortran runtime error: Actual string # length is shorter than the declared one for dummy argument 'which' (0/2)
For asymmetric matrices, we need to go to each end separately, so finding the domain bounds can be hard. Using row-standardisation has the nice feature of setting the upper bound to unity, and there are graph methods for finding out whether the lower bound is -1
.
======================================================
First we'll see how to get from sparse matrices to graphs. The mode of a symmetric matrix is "undirected"
by definition:
class(B) object.size(B) library(igraph) g1 <- graph.adjacency(B, mode="undirected") class(g1) object.size(g1)
We can also convert this graph pack to the same matrix, but note that get.adjacency
chooses a particular class of sparse matrix to be returned, so that the conversion process typically leads many matrices to fewer graph types, and back to fewer matrix types:
B1 <- get.adjacency(g1) class(B1) object.size(B1) all.equal(B, as(as(B1, "dgTMatrix"), "symmetricMatrix"))
A simple example of using igraph to do the same as an existing spdep function is Nicholas Lewin-Koh's n.comp.nb
from the early days of the package. It is useful to know whether an nb
object is divided up into separate subgraphs, and which entities are members of which such subgraph.
res <- n.comp.nb(col2) table(res$comp.id)
The same result can be obtained using the clusters
function in igraph:
c1 <- clusters(g1) c1$no == res$nc all.equal(c1$membership, res$comp.id) all.equal(c1$csize, c(table(res$comp.id)), check.attributes=FALSE)
The same holds for the row-standardised variant:
W <- as(nb2listw(col2, style="W", zero.policy=TRUE), "CsparseMatrix") g1W <- graph.adjacency(W, mode="directed", weighted="W") c1W <- clusters(g1W) all.equal(c1W$membership, res$comp.id)
Finding shortest paths between spatial entities across a given graph is a way to express closeness. If the graph is connected, that is that it is possible to traverse the graph edges from any node to any other, the longest shortest path is then a useful measure. In igraph, the is.connected
function tells us tells us that our graph is not connected, as we know from the figure above. The diameter measure is then the diameter of the largest component subgraph. Note that this generates an n
x n
matrix:
is.connected(g1) dg1 <- diameter(g1) dg1 sp_mat <- shortest.paths(g1) str(sp_mat)
If we do the same in spdep, using nblag
to a maximum number of lag orders - the diameter, but which is unknown in advance (the largest lag order for which the number of links is greater than zero), we run into the problem of how to represent missing neighbour information.
nbl10 <- nblag(col2, maxlag=10) vals <- sapply(nbl10, function(x) sum(card(x))) zero <- which(vals == 0) zero[1]-1
If we insert zero into the weights matrix where there is no connection using zero.policy=TRUE
, we generate a zero shortest path. If we are to create a matrix that matches the one produced by shortest.paths
, we need to set all these non-structural zeros to infinity (the length of the path between unconnected nodes), and re-instate structural zeros on the diagonal:
lmat <- lapply(nbl10[1:(zero[1]-1)], nb2mat, style="B", zero.policy=TRUE) mat <- matrix(0, n, n) for (i in seq(along=lmat)) mat = mat + i*lmat[[i]] mat[mat==0] <- Inf diag(mat) <- 0 all.equal(mat, sp_mat, check.attributes=FALSE)
Another area in which a graph representation might prove useful is in trying to establish the domain of the spatial coefficient when spatial weights are row-standardised. In that case by construction we know that the maximum eigenvalue is 1. If there are multiple blocks, that is graph components, where the numbers of nodes per block are greater than 1, then each will have a maximum eigenvalue of 1. The remaining problems are the numbers of zero eigenvalues (at least the singleton graph components), and whether any non-singleton component fulfills the condition termed by Smirnov and Anselin (2009) a cyclical matrix, for which the minimum eigenvalue is -1. The term cyclical appears to be used in many different ways, and it is not clear that its use here after Smirnov and Anselin (2009, pp. 2984-2985) indicates which meaning should be used to find the relevant graph function. The definition used here is that a block matrix (subgraph) is cyclical if: "for every location, every pair of its neighbours are not connected." That is, if w[i,j] and w[i,k] are greater than zero, w[j,k] must be zero to meet the condition.
The internal function find_q1_q2 returns the number of non-singleton components, and the number of these that meet this condition. It does this for each block/subgraph by testing the condition until it meets w[j,k] > 0, at which point it breaks. Smirnov and Anselin (2009) state that rook neighbours on a regular grid meet the condition:
nb_r <- cell2nb(7, 7, type="rook") nb_rW <- nb2listw(nb_r, style="W") spdep:::find_q1_q2(nb_rW)
One block/graph component is found, and this one meets the cyclical matrix condition, as also shown by the domain:
1/range(Re(eigenw(similar.listw(nb_rW))))
This does not apply to the spatial weights we have been using above, with two non-singleton components, neither meeting the cyclical matrix condition:
spdep:::find_q1_q2(nb_W) 1/range(Re(eigenw(similar.listw(nb_W))))
By construction, all two-node connected graph components also meet the condition, as the eigenvalues sum to zero, and the maximum is unity, so the minimum must be -1.
Any scripts or data that you put into this service are public.
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.