interp | R Documentation |
This function implements bivariate interpolation for irregularly spaced input data. Piecewise linear (=barycentric interpolation), bilinear or bicubic spline interpolation according to Akimas method is applied.
interp(x, y = NULL, z, xo = seq(min(x), max(x), length = nx),
yo = seq(min(y), max(y), length = ny),
linear = (method == "linear"), extrap = FALSE,
duplicate = "error", dupfun = NULL,
nx = 40, ny = 40, input="points", output = "grid",
method = "linear", deltri = "shull", h=0,
kernel="gaussian", solver="QR", degree=3,
baryweight=TRUE, autodegree=FALSE, adtol=0.1,
smoothpde=FALSE, akimaweight=TRUE, nweight=25,
na.rm=FALSE)
x |
vector of |
y |
vector of If left as NULL indicates that |
z |
vector of Missing values are not accepted by default, see parameter
|
xo |
If If |
yo |
If If |
input |
text, possible values are This is used to distinguish between regular and irregular gridded input data. |
output |
text, possible values are If In the case of |
linear |
logical, only for backward compatibility with Please use the new |
method |
text, possible methods are
|
extrap |
logical, indicates if extrapolation outside the convex hull is intended, this will not work for piecewise linear interpolation! |
duplicate |
character string indicating how to handle duplicate data points. Possible values are
|
dupfun |
a function, applied to duplicate points if
|
nx |
dimension of output grid in x direction |
ny |
dimension of output grid in y direction |
deltri |
triangulation method used, this argument may later be moved
into a control set together with others related to the spline
interpolation! Possible values are |
h |
bandwidth for partial derivatives estimation, compare |
kernel |
kernel for partial derivatives estimation, compare |
solver |
solver used in partial derivatives estimation, compare |
degree |
degree of local polynomial used for partial derivatives
estimation, compare |
baryweight |
calculate three partial derivatives estimators and return a barycentric weighted average. This increases the accuracy of Akima splines but the runtime is multplied by 3! |
autodegree |
try to reduce |
adtol |
tolerance used for autodegree |
smoothpde |
Use an averaged version of partial derivatives
estimates, by default simple average of Currently disabled by default (FALSE), underlying code still a bit experimental. |
akimaweight |
apply Akima weighting scheme on partial derivatives estimations instead of simply averaging |
nweight |
size of search neighbourhood for weighting scheme, default: 25 |
na.rm |
remove points where z= |
a list with 3 components:
x , y |
If If |
z |
If If If the input was a |
Please note that this function tries to be a replacement for the interp() function from the akima package. So it should be call compatible for most applications. It also offers additional tuning parameters, usually the default settings will fit. Please be aware that these additional parameters may change in the future as they are still under development.
Albrecht Gebhardt <albrecht.gebhardt@aau.at>, Roger Bivand <roger.bivand@nhh.no>
Moebius, A. F. (1827) Der barymetrische Calcul. Verlag v. Johann Ambrosius Barth, Leipzig, https://books.google.at/books?id=eFPluv_UqFEC&hl=de&pg=PR1#v=onepage&q&f=false
Franke, R., (1979). A critical comparison of some methods for interpolation of scattered data. Tech. Rep. NPS-53-79-003, Dept. of Mathematics, Naval Postgraduate School, Monterey, Calif.
Akima, H. (1978). A Method of Bivariate Interpolation and Smooth Surface Fitting for Irregularly Distributed Data Points. ACM Transactions on Mathematical Software 4, 148-164.
Akima, H. (1996). Algorithm 761: scattered-data surface fitting that has the accuracy of a cubic polynomial. ACM Transactions on Mathematical Software 22, 362–371.
interpp
### Use all datasets from Franke, 1979:
data(franke)
## x-y irregular grid points:
oldseed <- set.seed(42)
ni <- 64
xi <- runif(ni,0,1)
yi <- runif(ni,0,1)
xyi <- cbind(xi,yi)
## linear interpolation
fi <- franke.fn(xi,yi,1)
IL <- interp(xi,yi,fi,nx=80,ny=80,method="linear")
## prepare breaks and colors that match for image and contour:
breaks <- pretty(seq(min(IL$z,na.rm=TRUE),max(IL$z,na.rm=TRUE),length=11))
db <- breaks[2]-breaks[1]
nb <- length(breaks)
breaks <- c(breaks[1]-db,breaks,breaks[nb]+db)
colors <- terrain.colors(length(breaks)-1)
image(IL,breaks=breaks,col=colors,main="Franke function 1",
sub=paste("linear interpolation, ", ni,"points"))
contour(IL,add=TRUE,levels=breaks)
points(xi,yi)
## spline interpolation
fi <- franke.fn(xi,yi,1)
IS <- interp(xi,yi,fi,method="akima",
kernel="gaussian",solver="QR")
## prepare breaks and colors that match for image and contour:
breaks <- pretty(seq(min(IS$z,na.rm=TRUE),max(IS$z,na.rm=TRUE),length=11))
db <- breaks[2]-breaks[1]
nb <- length(breaks)
breaks <- c(breaks[1]-db,breaks,breaks[nb]+db)
colors <- terrain.colors(length(breaks)-1)
image(IS,breaks=breaks,col=colors,main="Franke function 1",
sub=paste("spline interpolation, ", ni,"points"))
contour(IS,add=TRUE,levels=breaks)
points(xi,yi)
## regular grid:
nx <- 8; ny <- 8
xg<-seq(0,1,length=nx)
yg<-seq(0,1,length=ny)
xx <- t(matrix(rep(xg,ny),nx,ny))
yy <- matrix(rep(yg,nx),ny,nx)
xyg<-expand.grid(xg,yg)
## linear interpolation
fg <- outer(xg,yg,function(x,y)franke.fn(x,y,1))
IL <- interp(xg,yg,fg,input="grid",method="linear")
## prepare breaks and colors that match for image and contour:
breaks <- pretty(seq(min(IL$z,na.rm=TRUE),max(IL$z,na.rm=TRUE),length=11))
db <- breaks[2]-breaks[1]
nb <- length(breaks)
breaks <- c(breaks[1]-db,breaks,breaks[nb]+db)
colors <- terrain.colors(length(breaks)-1)
image(IL,breaks=breaks,col=colors,main="Franke function 1",
sub=paste("linear interpolation, ", nx,"x",ny,"points"))
contour(IL,add=TRUE,levels=breaks)
points(xx,yy)
## spline interpolation
fg <- outer(xg,yg,function(x,y)franke.fn(x,y,1))
IS <- interp(xg,yg,fg,input="grid",method="akima",
kernel="gaussian",solver="QR")
## prepare breaks and colors that match for image and contour:
breaks <- pretty(seq(min(IS$z,na.rm=TRUE),max(IS$z,na.rm=TRUE),length=11))
db <- breaks[2]-breaks[1]
nb <- length(breaks)
breaks <- c(breaks[1]-db,breaks,breaks[nb]+db)
colors <- terrain.colors(length(breaks)-1)
image(IS,breaks=breaks,col=colors,main="Franke function 1",
sub=paste("spline interpolation, ", nx,"x",ny,"points"))
contour(IS,add=TRUE,levels=breaks)
points(xx,yy)
## apply interp to sp data:
require(sp)
## convert Akima data set to a sp object
data(akima)
asp <- SpatialPointsDataFrame(list(x=akima$x,y=akima$y),
data = data.frame(z=akima$z))
spplot(asp,"z")
## linear interpolation
spli <- interp(asp, z="z", method="linear")
## the result is again a SpatialPointsDataFrame:
spplot(spli,"z")
## now with spline interpolation, slightly higher resolution
spsi <- interp(asp, z="z", method="akima", nx=120, ny=120)
spplot(spsi,"z")
## now sp grids: reuse stuff from above
spgr <- SpatialPixelsDataFrame(list(x=c(xx),y=c(yy)),
data=data.frame(z=c(fg)))
spplot(spgr)
## linear interpolation
spli <- interp(spgr, z="z", method="linear", input="grid")
## the result is again a SpatialPointsDataFrame:
spplot(spli,"z")
## now with spline interpolation, slightly higher resolution
spsi <- interp(spgr, z="z", method="akima", nx=240, ny=240)
spplot(spsi,"z")
set.seed(oldseed)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.