SSJ  3.3.1
Stochastic Simulation in Java
Public Member Functions | Package Attributes | List of all members
SmoothingCubicSpline Class Reference

Represents a cubic spline with nodes at \((x_i, y_i)\) computed with the smoothing cubic spline algorithm of Schoenberg [42], [204] . More...

Inheritance diagram for SmoothingCubicSpline:
[legend]
Collaboration diagram for SmoothingCubicSpline:
[legend]

Public Member Functions

 SmoothingCubicSpline (double[] x, double[] y, double[] w, double rho)
 Constructs a spline with nodes at \((x_i, y_i)\), with weights \(w_i\) and smoothing factor \(\rho\) = rho. More...
 
 SmoothingCubicSpline (double[] x, double[] y, double rho)
 Constructs a spline with nodes at \((x_i, y_i)\), with weights \(= 1\) and smoothing factor \(\rho\) = rho. More...
 
double evaluate (double z)
 Evaluates and returns the value of the spline at \(z\). More...
 
double integral (double a, double b)
 Evaluates and returns the value of the integral of the spline from \(a\) to \(b\). More...
 
double derivative (double z)
 Evaluates and returns the value of the first derivative of the spline at \(z\). More...
 
double derivative (double z, int n)
 Evaluates and returns the value of the n-th derivative of the spline at \(z\). More...
 
double [] getX ()
 Returns the \(x_i\) coordinates for this spline. More...
 
double [] getY ()
 Returns the \(y_i\) coordinates for this spline. More...
 
double [] getWeights ()
 Returns the weights of the points. More...
 
double getRho ()
 Returns the smoothing factor used to construct the spline. More...
 
Polynomial [] getSplinePolynomials ()
 Returns a table containing all fitting polynomials. More...
 
int getFitPolynomialIndex (double x)
 Returns the index of \(P\), the umontreal.ssj.functions.Polynomial instance used to evaluate \(x\), in an ArrayList table instance returned by getSplinePolynomials(). More...
 

Package Attributes

double [] y
 
double [] weight
 

Detailed Description

Represents a cubic spline with nodes at \((x_i, y_i)\) computed with the smoothing cubic spline algorithm of Schoenberg [42], [204] .

A smoothing cubic spline is made of \(n+1\) cubic polynomials. The \(i\)th polynomial of such a spline, for \(i=1,…,n-1\), is defined as \(S_i(x)\) while the complete spline is defined as

\[ S(x) = S_i(x), \qquad\mbox{ for }x\in[x_{i-1}, x_i]. \]

For \(x<x_0\) and \(x>x_{n-1}\), the spline is not precisely defined, but this class performs extrapolation by using \(S_0\) and \(S_n\) linear polynomials. The algorithm which calculates the smoothing spline is a generalization of the algorithm for an interpolating spline. \(S_i\) is linked to \(S_{i+1}\) at \(x_{i+1}\) and keeps continuity properties for first and second derivatives at this point, therefore \(S_i(x_{i+1})=S_{i+1}(x_{i+1})\), \(S’_i(x_{i+1})=S’_{i+1}(x_{i+1})\) and \(S"_i(x_{i+1})=S"_{i+1}(x_{i+1})\).

The spline is computed with a smoothing parameter \(\rho\in[0, 1]\) which represents its accuracy with respect to the initial \((x_i, y_i)\) nodes. The smoothing spline minimizes

\[ L = \rho\sum_{i=0}^{n-1}{w_i}\left({y_i - S_i(x_i)}\right)^2 + (1-\rho)\int_{x_0}^{x_{n-1}}\left(S"(x)\right)^2dx \]

In fact, by setting \(\rho= 1\), we obtain the interpolating spline; and we obtain a linear function by setting \(\rho= 0\). The weights \(w_i>0\), which default to 1, can be used to change the contribution of each point in the error term. A large value \(w_i\) will give a large weight to the \(i\)th 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;
int m = 40;
double[] Xp = new double[m+1]; // Xp, Yp are spline points
double[] Yp = new double[m+1];
double h = (X[n-1] - X[0]) / m; // step
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
}

Constructor & Destructor Documentation

◆ SmoothingCubicSpline() [1/2]

SmoothingCubicSpline ( double []  x,
double []  y,
double []  w,
double  rho 
)

Constructs a spline with nodes at \((x_i, y_i)\), with weights \(w_i\) and smoothing factor \(\rho\) = rho.

The \(x_i\) must be sorted in increasing order.

Parameters
xthe \(x_i\) coordinates.
ythe \(y_i\) coordinates.
wthe weight for each point, must be \(> 0\).
rhothe smoothing parameter
Exceptions
IllegalArgumentExceptionif x, y and z do not have the same length, if rho has wrong value, or if the spline cannot be calculated.

◆ SmoothingCubicSpline() [2/2]

SmoothingCubicSpline ( double []  x,
double []  y,
double  rho 
)

Constructs a spline with nodes at \((x_i, y_i)\), with weights \(= 1\) and smoothing factor \(\rho\) = rho.

The \(x_i\) must be sorted in increasing order.

Parameters
xthe \(x_i\) coordinates.
ythe \(y_i\) coordinates.
rhothe smoothing parameter
Exceptions
IllegalArgumentExceptionif x and y do not have the same length, if rho has wrong value, or if the spline cannot be calculated.

Member Function Documentation

◆ derivative() [1/2]

double derivative ( double  z)

Evaluates and returns the value of the first derivative of the spline at \(z\).

Parameters
zargument of the spline.
Returns
value of first derivative.

Implements MathFunctionWithFirstDerivative.

◆ derivative() [2/2]

double derivative ( double  z,
int  n 
)

Evaluates and returns the value of the n-th derivative of the spline at \(z\).

Parameters
zargument of the spline.
norder of the derivative.
Returns
value of n-th derivative.

Implements MathFunctionWithDerivative.

◆ evaluate()

double evaluate ( double  z)

Evaluates and returns the value of the spline at \(z\).

Parameters
zargument of the spline.
Returns
value of spline.

Implements MathFunction.

◆ getFitPolynomialIndex()

int getFitPolynomialIndex ( double  x)

Returns the index of \(P\), the umontreal.ssj.functions.Polynomial instance used to evaluate \(x\), in an ArrayList table instance returned by getSplinePolynomials().

This index \(k\) gives also the interval in table X which contains the value \(x\) (i.e. such that \(x_k < x \leq x_{k+1}\)).

Returns
Index of the polynomial check with x in the Polynomial list returned by getSplinePolynomials

◆ getRho()

double getRho ( )

Returns the smoothing factor used to construct the spline.

Returns
the smoothing factor.

◆ getSplinePolynomials()

Polynomial [] getSplinePolynomials ( )

Returns a table containing all fitting polynomials.

Returns
Table containing the fitting polynomials.

◆ getWeights()

double [] getWeights ( )

Returns the weights of the points.

Returns
the weights of the points.

◆ getX()

double [] getX ( )

Returns the \(x_i\) coordinates for this spline.

Returns
the \(x_i\) coordinates.

◆ getY()

double [] getY ( )

Returns the \(y_i\) coordinates for this spline.

Returns
the \(y_i\) coordinates.

◆ integral()

double integral ( double  a,
double  b 
)

Evaluates and returns the value of the integral of the spline from \(a\) to \(b\).

Parameters
alower limit of integral.
bupper limit of integral.
Returns
value of integral.

Implements MathFunctionWithIntegral.


The documentation for this class was generated from the following file: