Guide Function Fit
Guide Function Fit
Package functionfit
Function fit utilities
Version: September 29, 2015
This package provides basic facilities for curve fitting and interpolation with polynomials
as, for example, least square fit and spline interpolation.
CONTENTS 1
Contents
PolInterp . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . .
LeastSquares . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . .
BSpline . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . .
SmoothingCubicSpline . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . .
PolInterp
Represents a polynomial that interpolates through a set of points. More specifically,
let (x0 , y0 ), . . . , (xn , yn ) be a set of points and p(x) the constructed polynomial of degree n.
Then, for i = 0, . . . , n, p(xi ) = yi .
package umontreal.iro.lecuyer.functionfit;
public class PolInterp extends Polynomial implements Serializable
Constructors
public PolInterp (double[] x, double[] y)
Constructs a new polynomial interpolating through the given points (x[0], y[0]), ...,
(x[n], y[n]). This constructs a polynomial of degree n from n+1 points.
Methods
public static double[] getCoefficients (double[] x, double[] y)
Computes and returns the coefficients the polynomial interpolating through the given points
(x[0], y[0]), ..., (x[n], y[n]). This polynomial has degree n and there are n+1 coefficients.
LeastSquares
This class implements different linear regression models, using the least squares method
to estimate the regression coefficients. Given input data xij and response yi , one want to
find the coefficients j that minimize the residuals of the form (using matrix notation)
r = min kY Xk2 ,
y i 0
k
X
!2
j xij
j=1
Sometimes, one wants to use a basis of general functions j (t) with a minimization of
the form
!2
k
X
X
r = min
yi
j j (ti ) .
j=1
For example, we could have j (t) = ej t or some other functions. In that case, one has to
choose the points ti at which to compute the basis functions, and use a method below with
xij = j (ti ).
package umontreal.iro.lecuyer.functionfit;
public class LeastSquares
Methods
public static double[] calcCoefficients (double[] X, double[] Y)
Computes the regression coefficients using the least squares method. This is a simple linear
regression with 2 regression coefficients, and . The model is
y = + x.
Given the n data points (Xi , Yi ), i = 0, 1, . . . , (n 1), the method computes and returns the
array [, ].
LeastSquares 4
Given the n data points (Xi , Yi ), i = 0, 1, . . . , (n 1), the method computes and returns the
array [0 , 1 , . . . , k ]. Restriction: n > k.
k1
X
j xj .
j=0
There are n data points Yi , Xij , i = 0, 1, . . . , (n 1), and each Xi is a k-dimensional point.
Given the response Y[i] and the regressor variables X[i][j], i = 0, 1, . . . , (n 1), j =
0, 1, . . . , (k 1), the method computes and returns the array [0 , 1 , . . . , k1 ]. Restriction:
n > k.
BSpline
Represents a B-spline with control points at (Xi , Yi ). Let Pi = (Xi , Yi ), for i = 0, . . . , n
1, be a control point and let tj , for j = 0, . . . , m 1 be a knot. A B-spline [1] of degree
p = m n 1 is a parametric curve defined as
P(t) =
n1
X
i=0
Here,
ti+p+1 t
t ti
Ni,p1 (t) +
Ni+1,p1 (t)
ti+p ti
ti+p+1 ti+1
1 for ti t ti+1 ,
Ni,0 (t) =
0 elsewhere.
Ni,p (t) =
This class provides methods to evaluate P(t) = (X(t), Y (t)) at any value of t, for a Bspline of any degree p 1. Note that the evaluate method of this class can be slow, since it
uses a root finder to determine the value of t for which X(t ) = x before it computes Y (t ).
package umontreal.iro.lecuyer.functionfit;
public class BSpline implements MathFunction
Constructors
public BSpline (final double[] x, final double[] y, final int degree)
Constructs a new uniform B-spline of degree degree with control points at (x[i], y[i]).
The knots of the resulting B-spline are set uniformly from x[0] to x[n-1].
Methods
public double[] getX()
Returns the Xi coordinates for this spline.
BSpline 6
n
X
Yi Si (Xi ) 2
i=0
Wi
This method uses the uniformly spaced method for interpolating points with a B-spline curve
and a uniformed clamped knot vector, as described in http://www.cs.mtu.edu/~shene/
COURSES/cs3621/NOTES/.
SmoothingCubicSpline
Represents a cubic spline with nodes at (xi , yi ) computed with the smoothing cubic spline
algorithm of Schoenberg [1, 2]. A smoothing cubic spline is made of n + 1 cubic polynomials.
The ith polynomial of such a spline, for i = 1, . . . , n1, is defined as Si (x) while the complete
spline is defined as
S(x) = Si (x),
for x [xi1 , xi ].
For x < x0 and x > xn1 , the spline is not precisely defined, but this class performs extrapolation by using S0 and Sn linear polynomials. The algorithm which calculates the smoothing
spline is a generalization of the algorithm for an interpolating spline. Si is linked to Si+1 at
xi+1 and keeps continuity properties for first and second derivatives at this point, therefore
00
0
(xi+1 ).
(xi+1 ) and Si00 (xi+1 ) = Si+1
Si (xi+1 ) = Si+1 (xi+1 ), Si0 (xi+1 ) = Si+1
The spline is computed with a smoothing parameter [0, 1] which represents its accuracy with respect to the initial (xi , yi ) nodes. The smoothing spline minimizes
L=
n1
X
xn1
wi (yi Si (xi )) + (1 )
i=0
(S 00 (x)) dx
x0
In fact, by setting = 1, we obtain the interpolating spline; and we obtain a linear function
by setting = 0. The weights wi > 0, which default to 1, can be used to change the
contribution of each point in the error term. A large value wi will give a large weight to the
ith point, so the spline will pass closer to it. Here is a small example that uses smoothing
splines:
int n;
double[] X = new double[n];
double[] Y = new double[n];
// here, fill arrays X and Y with n data points (x_i, y_i)
// The points must be sorted with respect to x_i.
double rho = 0.1;
SmoothingCubicSpline fit = new SmoothingCubicSpline(X, Y, rho);
int m = 40;
double[] Xp = new double[m+1];
double[] Yp = new double[m+1];
double h = (X[n-1] - X[0]) / m;
for (int i = 0; i <= m; i++) {
double z = X[0] + i * h;
Xp[i] = z;
Yp[i] = fit.evaluate(z);
}
// evaluate spline at z
SmoothingCubicSpline 8
package umontreal.iro.lecuyer.functionfit;
import umontreal.iro.lecuyer.functions.*;
import umontreal.iro.lecuyer.functions.Polynomial;
public class SmoothingCubicSpline implements MathFunction,
MathFunctionWithFirstDerivative, MathFunctionWithDerivative,
MathFunctionWithIntegral
Constructors
public SmoothingCubicSpline (double[] x, double[] y, double[] w,
double rho)
Constructs a spline with nodes at (xi , yi ), with weights wi and smoothing factor = rho.
The xi must be sorted in increasing order.
Methods
public double evaluate (double z)
Evaluates and returns the value of the spline at z.
SmoothingCubicSpline 9
REFERENCES 10
References
[1] C. de Boor. A Practical Guide to Splines. Number 27 in Applied Mathematical Sciences
Series. Springer-Verlag, New York, 1978.
[2] D. S. G. Pollock. Smoothing with cubic splines. Technical report, University of London,
Queen Mary and Westfield College, London, 1993.