Description Usage Arguments Details Value Author(s) References See Also Examples
The eigenw function returns a numeric vector of eigenvalues of 
the weights matrix generated from the spatial weights object listw. 
The eigenvalues are used to speed the computation of the Jacobian in 
spatial model estimation:
\log(\det[I - ρ W]) = ∑_{i=1}^{n}\log(1 - ρ λ_i)
where W is the n by n spatial weights matrix, and lambda[i] are the eigenvalues of W.
| 1 2 3 | eigenw(listw, quiet=NULL)
griffith_sone(P, Q, type="rook")
subgraph_eigenw(nb, glist=NULL, style="W", zero.policy=NULL, quiet=NULL)
 | 
| listw | a  | 
| quiet | default NULL, use global !verbose option value; set to FALSE for short summary | 
| P | number of columns in the grid (number of units in a horizontal axis direction) | 
| Q | number of rows in the grid (number of units in a vertical axis direction.) | 
| type | “rook” or “queen” | 
| nb | an object of class  | 
| glist | list of general weights corresponding to neighbours | 
| style | 
 | 
| zero.policy | default NULL, use global option value; if FALSE stop with error for any empty neighbour sets, if TRUE permit the weights list to be formed with zero-length weights vectors | 
The eigenw function computes the eigenvalues of a single spatial weights object. The griffith_sone function function may be used, following Ord and Gasim (for references see Griffith and Sone (1995)), to calculate analytical eigenvalues for binary rook or queen contiguous neighbours where the data are arranged as a regular P times Q grid. The subgraph_eigenw function may be used when there are multiple graph components, of which the largest may be handled as a dense matrix. Here the eigenvalues are computed for each subgraph in turn, and catenated to reconstruct the complete set. The functions may be used to provide pre-computed eigenvalues for spatial regression functions.
a numeric or complex vector of eigenvalues of the weights matrix generated from the spatial weights object.
Roger Bivand Roger.Bivand@nhh.no
Cliff, A. D., Ord, J. K. 1981 Spatial processes, Pion, p. 155; Ord, J. K. 1975 Estimation methods for models of spatial interaction, Journal of the American Statistical Association, 70, 120-126.; Griffith, D. A. and Sone, A. (1995). Trade-offs associated with normalizing constant computational simplifications for estimating spatial statistical models. Journal of Statistical Computation and Simulation, 51, 165-183.
| 1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 | data(oldcol)
W.eig <- eigenw(nb2listw(COL.nb, style="W"))
1/range(W.eig)
S.eig <- eigenw(nb2listw(COL.nb, style="S"))
1/range(S.eig)
B.eig <- eigenw(nb2listw(COL.nb, style="B"))
1/range(B.eig)
# cases for intrinsically asymmetric weights
crds <- cbind(COL.OLD$X, COL.OLD$Y)
k3 <- knn2nb(knearneigh(crds, k=3))
is.symmetric.nb(k3)
k3eig <- eigenw(nb2listw(k3, style="W"))
is.complex(k3eig)
rho <- 0.5
Jc <- sum(log(1 - rho * k3eig))
# complex eigenvalue Jacobian
Jc
# subgraphs
nc <- n.comp.nb(k3)
nc$nc
table(nc$comp.id)
k3eigSG <- subgraph_eigenw(k3, style="W")
all.equal(sort(k3eig), k3eigSG)
W <- as(nb2listw(k3, style="W"), "CsparseMatrix")
I <- diag(length(k3))
Jl <- sum(log(abs(diag(slot(lu(I - rho * W), "U")))))
# LU Jacobian equals complex eigenvalue Jacobian
Jl
all.equal(Re(Jc), Jl)
# wrong value if only real part used
Jr <- sum(log(1 - rho * Re(k3eig)))
Jr
all.equal(Jr, Jl)
# construction of Jacobian from complex conjugate pairs (Jan Hauke)
Rev <- Re(k3eig)[which(Im(k3eig) == 0)]
# real eigenvalues
Cev <- k3eig[which(Im(k3eig) != 0)]
pCev <- Cev[Im(Cev) > 0]
# separate complex conjugate pairs
RpCev <- Re(pCev)
IpCev <- Im(pCev)
# reassemble Jacobian
Jc1 <- sum(log(1 - rho*Rev)) + sum(log((1 - rho * RpCev)^2 + (rho^2)*(IpCev^2)))
all.equal(Re(Jc), Jc1)
# impact of omitted complex part term in real part only Jacobian
Jc2 <- sum(log(1 - rho*Rev)) + sum(log((1 - rho * RpCev)^2))
all.equal(Jr, Jc2)
# trace of asymmetric (WW) and crossprod of complex eigenvalues for APLE
sum(diag(W %*% W))
crossprod(k3eig)
# analytical regular grid eigenvalues
rg <- cell2nb(ncol=7, nrow=7, type="rook")
rg_eig <- eigenw(nb2listw(rg, style="B"))
rg_GS <- griffith_sone(P=7, Q=7, type="rook")
all.equal(rg_eig, rg_GS)
if (require(igraph)) {
B <- as(nb2listw(rg, style="B"), "CsparseMatrix")
f <- function(x, extra=NULL) {as.vector(B %*% x)}
res1 <- arpack(f, sym=TRUE, options=list(n=nrow(B), nev=1, ncv=8, which="LA",
 maxiter=200))
resn <- arpack(f, sym=TRUE, options=list(n=nrow(B), nev=1, ncv=8, which="SA", maxiter=200))
print(c(resn$value, res1$value))
#res <- arpack(f, sym=TRUE, options=list(n=nrow(B), nev=2, ncv=8, which="BE",
# maxiter=200))
#print(res$value)
# At line 558 of file dsaup2.f: Fortran runtime error:
# Index '9' of dimension 1 of array 'bounds' above upper bound of 8
print(all.equal(range(Re(rg_eig)), c(resn$value, res1$value)))
lw <- nb2listw(rg, style="W")
rg_eig <- eigenw(similar.listw(lw))
print(range(Re(rg_eig)))
W  <- as(lw, "CsparseMatrix")
f <- function(x, extra=NULL) {as.vector(W %*% x)}
print(arpack(f, sym=FALSE, options=list(n=nrow(W), nev=1, ncv=8, which="LR",
 maxiter=200))$value)
print(arpack(f, sym=FALSE, options=list(n=nrow(W), nev=1, ncv=8, which="SR",
 maxiter=200))$value)
}
 | 
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.