| Title: | Performs interpolation in multidimensional settings |
|---|---|
| Description: | This package implements multivariate extensions of the interpolation algorithms from Shepard, Donald (1968): "A two-dimensional interpolation function for irregularly-spaced data". Proceedings of the 1968 ACM National Conference. pp. 517<96>524. doi:10.1145/800186.810616. |
| Authors: | Andreas Maendle [aut, cre] |
| Maintainer: | Andreas Maendle <[email protected]> |
| License: | GPL-3 |
| Version: | 0.0.0.9000 |
| Built: | 2026-05-18 09:41:33 UTC |
| Source: | https://github.com/amaendle/mvInterpolation |
Roxygen2 allows you to write documentation in comment blocks co-located with code.
The only function Otherwise refer to the vignettes to see how to format the documentation.
Andreas Mändle, [email protected]
Shepard, Donald (1968). "A two-dimensional interpolation function for irregularly-spaced data". Proceedings of the 1968 ACM National Conference. pp. 517–524. doi:10.1145/800186.810616.
This function is basically a helper function for mvInterpolation::ti_get and not intended for direct use.
cosangmx(data, xinp, distance = NULL)cosangmx(data, xinp, distance = NULL)
data |
A matrix of data points |
xinp |
Coordinates of the fixed point of the angles; usually the point where interpolation has to be computed |
distance |
Distances between xinp and data can be given as an optional argument, if they have been coputed before, such that a new computation is not necessary |
The cosini of the pairwise angles between two data points with fixed point xinp
coords <- matrix(runif(10*4)*10, ncol=4) cosangmx(coords,rep(4,4))coords <- matrix(runif(10*4)*10, ncol=4) cosangmx(coords,rep(4,4))
Interpolation weights for inverse distance weighting.
di_get(data, xinp, p = 2)di_get(data, xinp, p = 2)
data |
A matrix of the points for which the function is known |
xinp |
Coordinates, where interpolation has to be computed |
p |
additional parameter as in the Shepard paper |
The weights for the IDW interpolation
mydata <- matrix(runif(10*4)*10, ncol=4) mydata <- cbind(mydata,abs(apply(mydata,1,sum)-3),abs(apply(mydata,1,prod)-4)) di_get(mydata, c(4,4,4,4))mydata <- matrix(runif(10*4)*10, ncol=4) mydata <- cbind(mydata,abs(apply(mydata,1,sum)-3),abs(apply(mydata,1,prod)-4)) di_get(mydata, c(4,4,4,4))
This function determines a search radius such that on average 7 points lie within the ball with the radius.
getR(data, d = NULL)getR(data, d = NULL)
data |
A matrix of the points for which the function is known |
d |
Dimension of the coordinate space |
The radius as a numeric
coords <- matrix(runif(10*4)*10, ncol=4) getR(coords)coords <- matrix(runif(10*4)*10, ncol=4) getR(coords)
Interpolation weights for inverse distance weighting of the nearest neighbors (4 to 10 neighbors, depending on the distance)
si_get(data, xinp, rad = getR)si_get(data, xinp, rad = getR)
data |
A matrix of the points for which the function is known |
xinp |
Coordinates, where interpolation has to be computed |
rad |
Function which gives the initial radius for the observations considered in the interpolation |
The weights for the IDW interpolation based on the nearest neighbors
mydata <- matrix(runif(10*4)*10, ncol=4) mydata <- cbind(mydata,abs(apply(mydata,1,sum)-3),abs(apply(mydata,1,prod)-4)) si_get(mydata, c(4,4,4,4))mydata <- matrix(runif(10*4)*10, ncol=4) mydata <- cbind(mydata,abs(apply(mydata,1,sum)-3),abs(apply(mydata,1,prod)-4)) si_get(mydata, c(4,4,4,4))
The drift is needed for interpolation by angular distance weighting extended by slopes. Interpolation weights are based on the nearest neighbors (4 to 10 neighbors, depending on the distance).
slopezi_get(data, xinp, rad = getR, vparam = 0.1)slopezi_get(data, xinp, rad = getR, vparam = 0.1)
data |
A matrix of the points for which the function is known |
xinp |
Coordinates, where interpolation has to be computed |
rad |
Function which gives the initial radius for the observations considered in the interpolation |
vparam |
Scalar, parameter for the slopes |
The weights for the ADW interpolation considering slopes and based on the nearest neighbors
mydata <- matrix(runif(10*4)*10, ncol=4) mydata <- cbind(mydata,abs(apply(mydata,1,sum)-3),abs(apply(mydata,1,prod)-4)) slopezi_get(mydata, c(4,4,4,4))mydata <- matrix(runif(10*4)*10, ncol=4) mydata <- cbind(mydata,abs(apply(mydata,1,sum)-3),abs(apply(mydata,1,prod)-4)) slopezi_get(mydata, c(4,4,4,4))
Interpolation weights for angular distance weighting of the nearest neighbors (4 to 10 neighbors, depending on the distance)
ti_get(data, xinp, rad = getR)ti_get(data, xinp, rad = getR)
data |
A matrix of the points for which the function is known |
xinp |
Coordinates, where interpolation has to be computed |
rad |
Function which gives the initial radius for the observations considered in the interpolation |
The weights for the ADW interpolation based on the nearest neighbors
mydata <- matrix(runif(10*4)*10, ncol=4) mydata <- cbind(mydata,abs(apply(mydata,1,sum)-3),abs(apply(mydata,1,prod)-4)) ti_get(mydata, c(4,4,4,4))mydata <- matrix(runif(10*4)*10, ncol=4) mydata <- cbind(mydata,abs(apply(mydata,1,sum)-3),abs(apply(mydata,1,prod)-4)) ti_get(mydata, c(4,4,4,4))
The drift is needed for interpolation by angular distance weighting extended by slopes. Interpolation weights include either all observations or are based on the nearest neighbors (4 to 10 neighbors, depending on the distance).
wintpl(data, xinp, weights, drift = 0, p = NULL, vparam = NULL)wintpl(data, xinp, weights, drift = 0, p = NULL, vparam = NULL)
data |
A matrix of the points for which the function is known |
xinp |
Coordinates, where interpolation has to be computed |
weights |
Vector of the weights for the interpolation |
drift |
Matrix of the drift for the interpolation |
p |
Scalar, parameter for the weights, currently not used |
vparam |
Scalar, parameter for the slopes, currently not used |
The interpoliation works for a multidimensional input space and also for multi-responses.
The interpolated value(s) at xinp
mydata <- matrix(runif(10*4)*10, ncol=4) mydata <- cbind(mydata,abs(apply(mydata,1,sum)-3),abs(apply(mydata,1,prod)-4)) inp <- rep(4,4) #Inverse distance weighting (based on all observations) wintpl(data=mydata, xinp=inp, weights=di_get(data=mydata, xinp=inp,2)) #Inverse distance weighting (based on nearest neighbors) wintpl(data=mydata, xinp=inp, weights=si_get(data=mydata, xinp=inp)^2) #Angular distance weighting wintpl(data=mydata, xinp=inp, weights=ti_get(data=mydata, xinp=inp)) #Angular distance weighting under consideration of slopes wintpl(data=mydata, xinp=inp, weights=ti_get(data=mydata, xinp=inp), drift=slopezi_get(mydata,inp))mydata <- matrix(runif(10*4)*10, ncol=4) mydata <- cbind(mydata,abs(apply(mydata,1,sum)-3),abs(apply(mydata,1,prod)-4)) inp <- rep(4,4) #Inverse distance weighting (based on all observations) wintpl(data=mydata, xinp=inp, weights=di_get(data=mydata, xinp=inp,2)) #Inverse distance weighting (based on nearest neighbors) wintpl(data=mydata, xinp=inp, weights=si_get(data=mydata, xinp=inp)^2) #Angular distance weighting wintpl(data=mydata, xinp=inp, weights=ti_get(data=mydata, xinp=inp)) #Angular distance weighting under consideration of slopes wintpl(data=mydata, xinp=inp, weights=ti_get(data=mydata, xinp=inp), drift=slopezi_get(mydata,inp))