splinefun {stats}R Documentation

Interpolating Splines

Description

Perform cubic (or Hermite) spline interpolation of given data points, returning either a list of points obtained by the interpolation or a function performing the interpolation.

Usage

splinefun(x, y = NULL,
          method = c("fmm", "periodic", "natural", "hyman",
                     "monoH.FC", "clamped", "canonical", "Catmull-Rom"),
          ties = mean, vals)

spline(x, y = NULL, n = 3*length(x), method = "fmm",
       xmin = min(x), xmax = max(x), xout, ties = mean)

splinefunH(x, y, m)

Arguments

x, y

vectors giving the coordinates of the points to be interpolated. Alternatively a single plotting structure can be specified: see xy.coords.

y must be increasing or decreasing for method = "hyman".

m

(for splinefunH()): vector of slopes m_i at the points (x_i,y_i); these together determine the Hermite “spline” which is piecewise cubic, (typically only) once differentiable continuously.

method

specifies the type of spline to be used. Possible values are "fmm", "natural", "periodic", "hyman", and only for splinefun() the Hermite spline based "monoH.FC", "clamped", "canonical", and "Catmull-Rom". Can be abbreviated.

n

if xout is left unspecified, interpolation takes place at n equally spaced points spanning the interval [xmin, xmax].

xmin, xmax

left-hand and right-hand endpoint of the interpolation interval (when xout is unspecified).

xout

an optional set of values specifying where interpolation is to take place.

ties

handling of tied x values. The string "ordered" or a function (or the name of a function) taking a single vector argument and returning a single number or a length-2 list of both, see approx and its ‘Details’ section, and the example below.

vals

for methods "clamped", the slopes at the two ends x[1] and x[length(x)]; for "canonical" the factors with which the.

Details

The inputs can contain missing values which are deleted, so at least one complete (x, y) pair is required. If method = "fmm", the spline used is that of Forsythe, Malcolm and Moler (an exact cubic is fitted through the four points at each end of the data, and this is used to determine the end conditions). Natural splines are used when method = "natural", and periodic splines when method = "periodic".

The method "monoH.FC" computes a monotone Hermite spline according to the method of \IFritsch and \ICarlson. It does so by determining slopes such that the Hermite spline (splinefunH()), determined by (x_i,y_i,m_i), is monotone (increasing or decreasing) iff the data are.

Method "hyman" computes a monotone cubic spline using \IHyman filtering of an method = "fmm" fit for strictly monotonic inputs.

These interpolation splines can also be used for extrapolation, that is prediction at points outside the range of x. Extrapolation makes little sense for method = "fmm"; for natural splines it is linear using the slope of the interpolating curve at the nearest data point.

Note that Hermite interpolation splines (via splinefunH()) are a more general class of functions than the others. They have more degrees of freedom with arbitrary slopes, and e.g., the natural interpolation spline (method = "natural") is the special case where the slopes are the divided differences.

Value

spline returns a list containing components x and y which give the ordinates where interpolation took place and the interpolated values.

splinefun returns a function with formal arguments x and deriv, the latter defaulting to zero. This function can be used to evaluate the interpolating cubic spline (deriv = 0), or its derivatives (deriv = 1, 2, 3) at the points x, where the spline function interpolates the data points originally specified. It uses data stored in its environment when it was created, the details of which are subject to change.

Warning

The value returned by splinefun contains references to the code in the current version of R: it is not intended to be saved and loaded into a different R session. This is safer in R >= 3.0.0.

Author(s)

R Core Team.

Simon Wood for the original code for Hyman filtering.

References

Becker, R. A., Chambers, J. M. and Wilks, A. R. (1988). The New S Language. Wadsworth & Brooks/Cole.

Dougherty, R. L., Edelman, A. and Hyman, J. M. (1989) Positivity-, monotonicity-, or convexity-preserving cubic and quintic Hermite interpolation. Mathematics of Computation, 52, 471–494. doi:10.1090/S0025-5718-1989-0962209-1.

Forsythe, G. E., Malcolm, M. A. and Moler, C. B. (1977). Computer Methods for Mathematical Computations. Wiley.

Fritsch, F. N. and Carlson, R. E. (1980). Monotone piecewise cubic interpolation. SIAM Journal on Numerical Analysis, 17, 238–246. doi:10.1137/0717021.

Hyman, J. M. (1983). Accurate monotonicity preserving cubic interpolation. SIAM Journal on Scientific and Statistical Computing, 4, 645–654. doi:10.1137/0904045.

https://en.wikiversity.org/wiki/Cubic_Spline_Interpolation on boundary conditions in addition to the natural ones;

https://en.wikipedia.org/wiki/Cubic_Hermite_spline for background on splinefunH() related methods "clamped", "canonical" and "Catmull-Rom".

See Also

approx and approxfun for constant and linear interpolation.

Package splines, especially interpSpline and periodicSpline for interpolation splines. That package also generates spline bases that can be used for regression splines.

smooth.spline for smoothing splines.

Examples

require(graphics)

op <- par(mfrow = c(2,1), mgp = c(2,.8,0), mar = 0.1+c(3,3,3,1))
n <- 9
x <- 1:n
y <- rnorm(n)
plot(x, y, main = paste("spline[fun](.) through", n, "points"))
lines(spline(x, y))
lines(spline(x, y, n = 201), col = 2)

y <- (x-6)^2
plot(x, y, main = "spline(.) -- 3 methods")
lines(spline(x, y, n = 201), col = 2)
lines(spline(x, y, n = 201, method = "natural"), col = 3)
lines(spline(x, y, n = 201, method = "periodic"), col = 4)
legend(6, 25, c("fmm","natural","periodic"), col = 2:4, lty = 1)

y <- sin((x-0.5)*pi)
f <- splinefun(x, y)
ls(envir = environment(f))
splinecoef <- get("z", envir = environment(f))
curve(f(x), 1, 10, col = "green", lwd = 1.5)
points(splinecoef, col = "purple", cex = 2)
curve(f(x, deriv = 1), 1, 10, col = 2, lwd = 1.5)
curve(f(x, deriv = 2), 1, 10, col = 2, lwd = 1.5, n = 401)
curve(f(x, deriv = 3), 1, 10, col = 2, lwd = 1.5, n = 401)
par(op)

fTrue <- function(x) exp(-3*x) * cos(4*pi*x)
x <- (0:20)/20
y <- fTrue(x)
f1  <- splinefun(x, y)# "fmm"
fn  <- splinefun(x, y, "natural")
f3  <- splinefun(x, y, "clamped", vals=c(-2, 1))
fCR <- splinefun(x, y, "Catmull-Rom")
f5  <- splinefun(x, y, "canonical", vals = 0.8)
x. <- seq(-1, 21, length.out = 1001)/20
matplot(x., cbind(f1(x.), fn(x.), f3(x.), fCR(x.), f5(x.)), type="l"); points(x,y)
ctrue <- adjustcolor("tomato", 1/4)
lines(x., fTrue(x.), col=ctrue, lwd=5)
legend("topright", col=c(palette()[1:5], ctrue), lty=1:5, lwd=c(1,1,1,1,1, 5), bty="n",
       c("fmm", "natural", "clamped(-2,1)", "Catmull-Rom","canonical(0.8)", "true f()"))

## Derivatives :
curve(f1 (x, deriv = 1), n = 1001); abline(h = 0, lty=3)
curve(fn (x, deriv = 1), add=TRUE, col=adjustcolor(2, 1/2), lwd=2.5, n = 1001)
curve(f3 (x, deriv = 1), add=TRUE, col=adjustcolor(3, 1/2), lwd=2.5, n = 1001)
curve(fCR(x, deriv = 1), add=TRUE, col=adjustcolor(4, 1/2), lwd=2.5, n = 1001)
curve(f5 (x, deriv = 1), add=TRUE, col=adjustcolor(5, 1/2), lwd=1.5, n = 1001)
legend("topright", col=c("black", adjustcolor(2:5, 1/2)), bty="n", lwd=c(1, rep(2.5,4)),
       c("fmm", "natural", "clamped(-2,1)", "Catmull-Rom","canonical(0.8)"))


## 2nd derivatives:  piecewise linear, only first three are continuous
curve(f1 (x, deriv = 2), n = 1001); abline(h = 0, lty=3)
curve(fn (x, deriv = 2), add=TRUE, col=adjustcolor(2, 1/2), lwd=2.5, n = 1001)
curve(f3 (x, deriv = 2), add=TRUE, col=adjustcolor(3, 1/2), lwd=2.5, n = 1001)
curve(fCR(x, deriv = 2), add=TRUE, col=adjustcolor(4, 1/2), lwd=2.5, n = 1001)
legend("topright", col=c("black", adjustcolor(2:5, 1/2)), bty="n", lwd=c(1, rep(2.5,4)),
       c("fmm", "natural", "clamped(-2,1)", "Catmull-Rom"))


## Manual spline evaluation --- demo the coefficients :
.x <- splinecoef$x
u <- seq(3, 6, by = 0.25)
(ii <- findInterval(u, .x))
dx <- u - .x[ii]
f.u <- with(splinecoef,
            y[ii] + dx*(b[ii] + dx*(c[ii] + dx* d[ii])))
stopifnot(all.equal(f(u), f.u))

## An example with ties (non-unique  x values):
set.seed(1); x <- round(rnorm(30), 1); y <- sin(pi * x) + rnorm(30)/10
plot(x, y, main = "spline(x,y)  when x has ties")
lines(spline(x, y, n = 201), col = 2)
## visualizes the non-unique ones:
tx <- table(x); mx <- as.numeric(names(tx[tx > 1]))
ry <- matrix(unlist(tapply(y, match(x, mx), range, simplify = FALSE)),
             ncol = 2, byrow = TRUE)
segments(mx, ry[, 1], mx, ry[, 2], col = "blue", lwd = 2)

## Another example with sorted x, but ties:
set.seed(8); x <- sort(round(rnorm(30), 1)); y <- round(sin(pi * x) + rnorm(30)/10, 3)
summary(diff(x) == 0) # -> 7 duplicated x-values
str(spline(x, y, n = 201, ties="ordered")) # all '$y' entries are NaN
## The default (ties=mean) is ok, but most efficient to use instead is
sxyo <- spline(x, y, n = 201, ties= list("ordered", mean))
sapply(sxyo, summary)# all fine now
plot(x, y, main = "spline(x,y, ties=list(\"ordered\", mean))  for when x has ties")
lines(sxyo, col="blue")

## An example of monotone interpolation
n <- 20
set.seed(11)
x. <- sort(runif(n)) ; y. <- cumsum(abs(rnorm(n)))
plot(x., y.)
curve(splinefun(x., y.)(x), add = TRUE, col = 2, n = 1001)
curve(splinefun(x., y., method = "monoH.FC")(x), add = TRUE, col = 3, n = 1001)
curve(splinefun(x., y., method = "hyman")   (x), add = TRUE, col = 4, n = 1001)
legend("topleft",
       paste0("splinefun( \"", c("fmm", "monoH.FC", "hyman"), "\" )"),
       col = 2:4, lty = 1, bty = "n")

## and one from Fritsch and Carlson (1980), Dougherty et al (1989)
x. <- c(7.09, 8.09, 8.19, 8.7, 9.2, 10, 12, 15, 20)
f <- c(0, 2.76429e-5, 4.37498e-2, 0.169183, 0.469428, 0.943740,
       0.998636, 0.999919, 0.999994)
s0 <- splinefun(x., f)
sn <- splinefun(x., f, method = "natural")
s1 <- splinefun(x., f, method = "monoH.FC")
s2 <- splinefun(x., f, method = "hyman")
plot(x., f, ylim = c(-0.2, 1.2))
curve(s0(x), add = TRUE, col = 2, n = 1001) -> m0
curve(sn(x), add = TRUE, col = 3, n = 1001)
curve(s1(x), add = TRUE, col = 4, n = 1001)
curve(s2(x), add = TRUE, col = 5, n = 1001)
legend("right",
       paste0("splinefun( \"", c("fmm", "natural", "monoH.FC", "hyman"), "\" )"),
       col = 2:5, lty = 1, bty = "n")

## "fmm" has continuous (piecewise linear) 2nd derivative :
curve(s0(x, deriv=2), xlim = range(x.), n=1000); abline(h=0, lty=3)
## but the monotone (hermite) splines do *not*:
cl <- adjustcolor(1:3, 1/2)
x.. <- sort(outer(x., (-1:1)/1000, "+"))
matplot(x.., cbind(s0(x.., 2), s1(x.., 2), s2(x.., 2)),
        type = "o", pch=1:3, col = cl, lty=1, lwd=2, cex=3/4)
legend("topright",
       paste0("splinefun( \"", c("fmm", "monoH.FC", "hyman"), "\" )"),
       pch=1:3, col = cl, lty=1, lwd=2, pt.cex=3/4, bty = "n")

## they seem identical, but are not quite:
xx <- m0$x
plot(xx, s1(xx) - s2(xx), type = "l",  col = 2, lwd = 2,
     main = "Difference   monoH.FC - hyman"); abline(h = 0, lty = 3)

x <- xx[xx < 10.2] ## full range: x <- xx .. does not show enough
ccol <- adjustcolor(2:4, 0.8)
matplot(x, cbind(s0(x, deriv = 2), s1(x, deriv = 2), s2(x, deriv = 2))^2,
        lwd = 2, col = ccol, type = "l", ylab = quote({{f*second}(x)}^2),
        main = expression({{f*second}(x)}^2 ~" for the three 'splines'"))
legend("topright",
       paste0("splinefun( \"", c("fmm", "monoH.FC", "hyman"), "\" )"),
       lwd = 2, col  =  ccol, lty = 1:3, bty = "n")
## --> "hyman" has slightly smaller  Integral f''(x)^2 dx  than "FC",
## here, and both are 'much worse' than the regular fmm spline.

[Package stats version 4.4.0 Index]