smbsp {bspline} | R Documentation |
Smoothing B-spline of order n >= 0
Description
Optimize smoothing B-spline coefficients (smbsp) and knot positions (fitsmbsp) such that residual squared sum is minimized for all y columns.
Usage
smbsp(
x,
y,
n = 3L,
xki = NULL,
nki = 1L,
lieq = NULL,
monotone = 0,
positive = 0,
mat = NULL,
estSD = FALSE,
tol = 1e-10
)
fitsmbsp(
x,
y,
n = 3L,
xki = NULL,
nki = 1L,
lieq = NULL,
monotone = 0,
positive = 0,
control = list(),
estSD = FALSE,
tol = 1e-10
)
Arguments
x |
Numeric vector, abscissa points |
y |
Numeric vector or matrix or data.frame, ordinate values to be smoothed (one set per column in case of matrix or data.frame) |
n |
Integer scalar, polynomial order of B-splines (3 by default) |
xki |
Numeric vector, strictly internal B-spline knots, i.e. lying strictly
inside of |
nki |
Integer scalar, internal knot number (1 by default). When
nki==0, it corresponds to polynomial regression. If |
lieq |
List, equality constraints to respect by the smoothing spline,
one list item per y column. By default (NULL), no constraint is imposed.
Constraints are given as a 2-column matrix |
monotone |
Numeric scalar or vector, if |
positive |
Numeric scalar, if |
mat |
Numeric matrix of basis vectors, if NULL it is recalculated
by |
estSD |
Logical scalar, if TRUE, indicates to calculate: SD of each y column, covariance matrix and SD of spline coefficients. All these values can be retrieved with bsppar() call (FALSE by default). These estimations are made under assumption that all y points have uncorrelated noise. Optional constraints are not taken into account of SD. |
tol |
Numerical scalar, relative tolerance for small singular values
that should be considered as 0 if |
control |
List, passed through to |
Details
If constraints are set, we use nlsic::lsie_ln()
to solve a
least squares
problem with equality constraints in least norm sens for each y column.
Otherwise, nlsic::ls_ln_svd()
is used for the whole y matrix.
The solution of least squares problem is a vector of B-splines coefficients qw
,
one vector per y
column. These vectors are used to define B-spline function
which is returned as the result.
NB. When nki >= length(x)-n-1
(be it from direct setting or calculated
from length(xki)
), it corresponds
to spline interpolation, i.e. the resulting spline will pass
exactly by (x,y) points (well, up to numerical precision).
Border and external knots are fixed, only strictly internal knots can move during optimization. The optimization process is constrained to respect a minimal distance between knots as well as to bound them to x range. This is done to avoid knots getting unsorted during iterations and/or going outside of a meaningful range.
Value
Function, smoothing B-splines
respecting optional constraints (generated by par2bsp()
).
See Also
bsppar
for retrieving parameters of B-spline functions; par2bsp
for generating B-spline function; iknots
for estimation of knot positions
Examples
x=seq(0, 1, length.out=11)
y=sin(pi*x)+rnorm(x, sd=0.1)
# constraint B-spline to be 0 at the interval ends
fsm=smbsp(x, y, nki=1, lieq=list(rbind(c(0, 0), c(1, 0))))
# check parameters of found B-splines
bsppar(fsm)
plot(x, y) # original "measurements"
# fine grained x
xfine=seq(0, 1, length.out=101)
lines(xfine, fsm(xfine)) # fitted B-splines
lines(xfine, sin(pi*xfine), col="blue") # original function
# visualize knot positions
xk=bsppar(fsm)$xk
points(xk, fsm(xk), pch="x", col="red")
# fit broken line with linear B-splines
x1=seq(0, 1, length.out=11)
x2=seq(1, 3, length.out=21)
x3=seq(3, 4, length.out=11)
y1=x1+rnorm(x1, sd=0.1)
y2=-2+3*x2+rnorm(x2, sd=0.1)
y3=4+x3+rnorm(x3, sd=0.1)
x=c(x1, x2, x3)
y=c(y1, y2, y3)
plot(x, y)
f=fitsmbsp(x, y, n=1, nki=2)
lines(x, f(x))