gwr.basic | R Documentation |
This function implements basic GWR
gwr.basic(formula, data, regression.points, bw, kernel="bisquare",
adaptive=FALSE, p=2, theta=0, longlat=F,dMat,F123.test=F,cv=F, W.vect=NULL,
parallel.method=FALSE,parallel.arg=NULL)
## S3 method for class 'gwrm'
print(x, ...)
formula |
Regression model formula of a formula object |
data |
a Spatial*DataFrame, i.e. SpatialPointsDataFrame or SpatialPolygonsDataFrame as defined in package sp, or a sf object defined in package sf |
regression.points |
a Spatial*DataFrame object, i.e. SpatialPointsDataFrame or SpatialPolygonsDataFrame as defined in package sp; Note that no diagnostic information will returned if it is assigned |
bw |
bandwidth used in the weighting function, possibly calculated by bw.gwr;fixed (distance) or adaptive bandwidth(number of nearest neighbours) |
kernel |
function chosen as follows: gaussian: wgt = exp(-.5*(vdist/bw)^2); exponential: wgt = exp(-vdist/bw); bisquare: wgt = (1-(vdist/bw)^2)^2 if vdist < bw, wgt=0 otherwise; tricube: wgt = (1-(vdist/bw)^3)^3 if vdist < bw, wgt=0 otherwise; boxcar: wgt=1 if dist < bw, wgt=0 otherwise |
adaptive |
if TRUE calculate an adaptive kernel where the bandwidth (bw) corresponds to the number of nearest neighbours (i.e. adaptive distance); default is FALSE, where a fixed kernel is found (bandwidth is a fixed distance) |
p |
the power of the Minkowski distance, default is 2, i.e. the Euclidean distance |
theta |
an angle in radians to rotate the coordinate system, default is 0 |
longlat |
if TRUE, great circle distances will be calculated |
dMat |
a pre-specified distance matrix, it can be calculated by the function |
F123.test |
If TRUE, conduct three seperate F-tests according to Leung et al. (2000). |
cv |
if TRUE, cross-validation data will be calculated and returned in the output Spatial*DataFrame |
W.vect |
default NULL, if given it will be used to weight the distance weighting matrix |
x |
an object of class “gwrm”, returned by the function |
parallel.method |
FALSE as default, and the calibration will be conducted traditionally via the serial technique, "omp": multi-thread technique with the OpenMP API, "cluster": multi-process technique with the parallel package, "cuda": parallel computing technique with CUDA |
parallel.arg |
if parallel.method is not FALSE, then set the argument by following: if parallel.method is "omp", parallel.arg refers to the number of threads used, and its default value is the number of cores - 1; if parallel.method is "cluster", parallel.arg refers to the number of R sessions used, and its default value is the number of cores - 1; if parallel.method is "cuda", parallel.arg refers to the number of calibrations included in each group, but note a too large value may cause the overflow of GPU memory. |
... |
arguments passed through (unused) |
A list of class “gwrm”:
GW.arguments |
a list class object including the model fitting parameters for generating the report file |
GW.diagnostic |
a list class object including the diagnostic information of the model fitting |
lm |
an object of class inheriting from “lm”, see lm. |
SDF |
a SpatialPointsDataFrame (may be gridded), or SpatialPolygonsDataFrame object (see package “sp”), or sf object (see package “sf”) integrated with regression.points, GWR coefficient estimates, y value,predicted values, coefficient standard errors and t-values in its "data" slot. |
timings |
starting and ending time. |
this.call |
the function call used. |
Ftest.res |
results of Leung's F tests when F123.test is TRUE. |
Requirements of using CUDA for high-performence computation in GWR functions:
To run GWR-CUDA (i.e. parallel.method is pecified as “cuda”) with
gwr.basic
, bw.gwr
and gwr.model.selection
,
the following conditions are required:
1. There is at least one NVIDIA GPU supporting CUDA equipped on user's computer.
2. CUDA (>10.2) are installed with the environment variable 'CUDA_HOME' set properly.
3. The package should re-built from source. - For Linux user, 'GWmodelCUDA' will be automatically built if CUDA toolkit could be detected by the complier. - For Windows user, 'GWmodelCUDA.dll' and 'GWmodelCUDA.lib' will be automatically downloaded; however, we would recommend users to build the 'GWmodelCUDA' library manually to avoid some potentially unknown issues, and an 'CMakeLists.txt' file is provided for this procedure.
If any condition above is not satisfied, the GWR-CUDA will not work even though the “parallel.method” is specified as “cuda”.
Binbin Lu binbinlu@whu.edu.cn
Brunsdon, C, Fotheringham, S, Charlton, M (1996), Geographically Weighted Regression: A Method for Exploring Spatial Nonstationarity. Geographical Analysis 28(4):281-298
Charlton, M, Fotheringham, S, and Brunsdon, C (2007), GWR3.0, http://gwr.nuim.ie/.
Fotheringham S, Brunsdon, C, and Charlton, M (2002), Geographically Weighted Regression: The Analysis of Spatially Varying Relationships, Chichester: Wiley.
Leung, Y, Mei, CL, and Zhang, WX (2000), Statistical tests for spatial nonstationarity based on the geographically weighted regression model. Environment and Planning A, 32, 9-32.
Lu, B, Charlton, M, Harris, P, Fotheringham, AS (2014) Geographically weighted regression with a non-Euclidean distance metric: a case study using hedonic house price data. International Journal of Geographical Information Science 28(4): 660-681
OpenMP: https://www.openmp.org/
CUDA: https://developer.nvidia.com/cuda-zone
R Core Team (2020). R: A language and environment for statistical computing. R Foundation for Statistical Computing, Vienna, Austria. https://www.R-project.org/.
data(LondonHP)
DM<-gw.dist(dp.locat=coordinates(londonhp))
##Compare the time consumed with and without a specified distance matrix
## Not run:
system.time(gwr.res<-gwr.basic(PURCHASE~FLOORSZ, data=londonhp, bw=1000,
kernel = "gaussian"))
system.time(DM<-gw.dist(dp.locat=coordinates(londonhp)))
system.time(gwr.res<-gwr.basic(PURCHASE~FLOORSZ, data=londonhp, bw=1000,
kernel = "gaussian", dMat=DM))
## specify an optimum bandwidth by cross-validation appraoch
bw1<-bw.gwr(PURCHASE~FLOORSZ, data=londonhp, kernel = "gaussian",dMat=DM)
gwr.res1<-gwr.basic(PURCHASE~FLOORSZ, data=londonhp, bw=bw1,kernel = "gaussian",
dMat=DM)
gwr.res1
## End(Not run)
data(LondonBorough)
nsa = list("SpatialPolygonsRescale", layout.north.arrow(), offset = c(561900,200900),
scale = 500, col=1)
## Not run:
if(require("RColorBrewer"))
{
mypalette<-brewer.pal(6,"Spectral")
x11()
spplot(gwr.res1$SDF, "FLOORSZ", key.space = "right", cex=1.5, cuts=10,
ylim=c(155840.8,200933.9), xlim=c(503568.2,561957.5),
main="GWR estimated coefficients for FLOORSZ with a fixed bandwidth",
col.regions=mypalette, sp.layout=list(nsa, londonborough))}
## End(Not run)
## Not run:
bw2<-bw.gwr(PURCHASE~FLOORSZ,approach="aic",adaptive=TRUE, data=londonhp,
kernel = "gaussian", dMat=DM)
gwr.res2<-gwr.basic(PURCHASE~FLOORSZ, data=londonhp, bw=bw2,adaptive=TRUE,
kernel = "gaussian", dMat=DM)
gwr.res2
if(require("RColorBrewer"))
{
x11()
spplot(gwr.res2$SDF, "FLOORSZ", key.space = "right", cex=1.5, cuts=10,
ylim=c(155840.8,200933.9), xlim=c(503568.2,561957.5),
main="GWR estimated coefficients for FLOORSZ with an adaptive bandwidth",
col.regions=mypalette, sp.layout=list(nsa,londonborough))}
## End(Not run)
## Not run:
############HP-GWR test code
simulate.data.generator <- function(data.length) {
x1 <- rnorm(data.length)
x2 <- rnorm(data.length)
x3 <- rnorm(data.length)
lon <- rnorm(data.length, mean = 533200, sd = 10000)
lat <- rnorm(data.length, mean = 159400, sd = 10000)
y <- x1 + 5 * x2 + 2.5 * x3 + rnorm(data.length)
simulate.data <- data.frame(y = y, x1 = x1, x2 = x2, x3 = x3, lon = lon, lat = lat)
coordinates(simulate.data) <- ~ lon + lat
names(simulate.data)
return(simulate.data)
}
simulate.data <- simulate.data.generator(10000)
adaptive = TRUE
## GWR (not parallelized)
bw.CV.s <- bw.gwr(data = simulate.data, formula = y ~ x1 + x2 + x3, approach="CV",
kernel = "gaussian", adaptive = adaptive, parallel.method = FALSE)
model.s <- gwr.model.selection(DeVar = "y", InDeVars = c("x1", "x2", "x3"), data = simulate.data,
bw = bw.CV.s, approach="AIC", kernel = "gaussian", adaptive = T,
parallel.method = FALSE)
system.time(
betas.s <- gwr.basic(data = simulate.data, formula = y ~ x1 + x2 + x3, bw = bw.CV.s,
kernel = "gaussian", adaptive = TRUE)
)
## GWR-Omp
bw.CV.omp <- bw.gwr(data = simulate.data, formula = y ~ x1 + x2 + x3, approach="CV",
kernel = "gaussian", adaptive = adaptive, parallel.method = "omp")
model.omp <- gwr.model.selection(DeVar = "y", InDeVars = c("x1", "x2", "x3"), data = simulate.data,
bw = bw.CV.omp, approach="AIC", kernel = "gaussian", adaptive = T,
parallel.method = "omp")
system.time(
betas.omp <- gwr.basic(data = simulate.data, formula = y ~ x1 + x2 + x3, bw = bw.CV.omp,
kernel = "gaussian", adaptive = T, parallel.method = "omp"))
## GWR-CUDA
bw.CV.cuda <- bw.gwr(data = simulate.data, formula = y ~ x1 + x2 + x3, approach="CV",
kernel = "gaussian", adaptive = adaptive, parallel.method = "cuda",
parallel.arg = 6*16)
model.cuda <- gwr.model.selection(DeVar = "y", InDeVars = c("x1", "x2", "x3"),
data = simulate.data, bw = bw.CV.cuda, approach="AIC",
kernel = "gaussian", adaptive = T,
parallel.method = "cuda", parallel.arg = 6*16)
system.time(
betas.cuda <- gwr.basic(data = simulate.data, formula = y ~ x1 + x2 + x3, bw = bw.CV.cuda,
kernel = "gaussian", adaptive = T, parallel.method = "cuda",
parallel.arg = 6*8))
## End(Not run)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.