Public API
API Summary
| Function | Description |
|---|---|
prepare | Prepare the spline by constructing and factoring a Gram matrix of the interpolation problem. |
construct | Construct the spline by calculating its coefficients. |
interpolate | Prepare and construct the spline. |
evaluate | Evaluate the spline value at the required locations |
evaluate_one | Evaluate the spline value at the required location |
evaluate_gradient | Evaluate gradient of the spline at the required location. |
evaluate_derivative | Evaluate the 1D spline derivative at the required location. |
estimate_accuracy | Estimate accuracy of the function interpolation result. |
estimate_cond | Estimate the Gram matrix 1-norm condition number. |
estimate_epsilon | Estimate the 'scaling parameter' of Bessel potential space the spline being built in. |
get_epsilon | Get the 'scaling parameter' of Bessel potential space the normal spline was built in. |
get_cond | Get the Gram matrix spectral condition number. |
Functions
NormalSmoothingSplines.prepare — Functionprepare(nodes::Matrix{T}, kernel::RK = RK_H0()) where {T <: AbstractFloat, RK <: ReproducingKernel_0}
Prepare the spline by constructing and factoring a Gram matrix of the interpolation problem. Initialize the NormalSpline object.
Arguments
nodes: The function value nodes. This should be ann×n_1matrix, wherenis dimension of the sampled space andn_1is the number of function value nodes. It means that each column in the matrix defines one node.kernel: reproducing kernel of Bessel potential space the normal spline is constructed in. It must be a struct object of the following type:RK_H0if the spline is constructing as a continuous function,RK_H1if the spline is constructing as a differentiable function,RK_H2if the spline is constructing as a twice differentiable function.
Return: the partly initialized NormalSpline object that must be passed to construct function in order to complete the spline initialization.
prepare(nodes::Matrix{T}, d_nodes::Matrix{T}, es::Matrix{T}, kernel::RK = RK_H1()) where {T <: AbstractFloat, RK <: ReproducingKernel_1}
Prepare the spline by constructing and factoring a Gram matrix of the interpolation problem. Initialize the NormalSpline object.
Arguments
nodes: The function value nodes. This should be ann×n_1matrix, wherenis dimension of the sampled space andn_1is the number of function value nodes. It means that each column in the matrix defines one node.d_nodes: The function directional derivatives nodes. This should be ann×n_2matrix, wherenis dimension of the sampled space andn_2is the number of function directional derivative nodes.es: Directions of the function directional derivatives. This should be ann×n_2matrix, wherenis dimension of the sampled space andn_2is the number of function directional derivative nodes. It means that each column in the matrix defines one direction of the function directional derivative.kernel: reproducing kernel of Bessel potential space the normal spline is constructed in. It must be a struct object of the following type:RK_H1if the spline is constructing as a differentiable function,RK_H2if the spline is constructing as a twice differentiable function.
Return: the partly initialized NormalSpline object that must be passed to construct function in order to complete the spline initialization.
prepare(nodes::Vector{T}, kernel::RK = RK_H0()) where {T <: AbstractFloat, RK <: ReproducingKernel_0}
Prepare the 1D spline by constructing and factoring a Gram matrix of the interpolation problem. Initialize the NormalSpline object.
Arguments
nodes: The function value nodes.kernel: reproducing kernel of Bessel potential space the normal spline is constructed in. It must be a struct object of the following type:RK_H0if the spline is constructing as a continuous function,RK_H1if the spline is constructing as a differentiable function,RK_H2if the spline is constructing as a twice differentiable function.
Return: the partly initialized NormalSpline object that must be passed to construct function in order to complete the spline initialization.
prepare(nodes::Vector{T}, d_nodes::Vector{T}, kernel::RK = RK_H1()) where {T <: AbstractFloat, RK <: ReproducingKernel_1}
Prepare the 1D normal spline by constructing and factoring a Gram matrix of the interpolation problem. Initialize the NormalSpline object.
Arguments
nodes: The function value nodes.d_nodes: The function derivatives nodes.kernel: reproducing kernel of Bessel potential space the normal spline is constructed in. It must be a struct object of the following type:RK_H1if the spline is constructing as a differentiable function,RK_H2if the spline is constructing as a twice differentiable function.
Return: the partly initialized NormalSpline object that must be passed to construct function in order to complete the spline initialization.
NormalSmoothingSplines.construct — Functionconstruct(spline::NormalSpline{T, RK}, values::Vector{T}) where {T <: AbstractFloat, RK <: ReproducingKernel_0}
Construct the spline by calculating its coefficients and completely initializing the NormalSpline object.
Arguments
spline: the partly initializedNormalSplineobject returned bypreparefunction.values: function values atnodesnodes.
Return: the completely initialized NormalSpline object that can be passed to evaluate function to interpolate the data to required points.
construct(spline::NormalSpline{T, RK}, values::Vector{T}, d_values::Vector{T}) where {T <: AbstractFloat, RK <: ReproducingKernel_1}
Construct the spline by calculating its coefficients and completely initializing the NormalSpline object.
Arguments
spline: the partly initializedNormalSplineobject returned bypreparefunction.values: function values atnodesnodes.d_values: function directional derivative values atd_nodesnodes.
Return: the completely initialized NormalSpline object that can be passed to evaluate function to interpolate the data to required points.
NormalSmoothingSplines.interpolate — Functioninterpolate(nodes::Matrix{T}, values::Vector{T}, kernel::RK = RK_H0()) where {T <: AbstractFloat, RK <: ReproducingKernel_0}
Prepare and construct the spline.
Arguments
nodes: The function value nodes. This should be ann×n_1matrix, wherenis dimension of the sampled space andn_1is the number of function value nodes. It means that each column in the matrix defines one node.values: function values atnodesnodes.kernel: reproducing kernel of Bessel potential space the normal spline is constructed in. It must be a struct object of the following type:RK_H0if the spline is constructing as a continuous function,RK_H1if the spline is constructing as a differentiable function,RK_H2if the spline is constructing as a twice differentiable function.
Return: the completely initialized NormalSpline object that can be passed to evaluate function.
interpolate(nodes::Matrix{T}, values::Vector{T}, d_nodes::Matrix{T}, es::Matrix{T}, d_values::Vector{T}, kernel::RK = RK_H1()) where {T <: AbstractFloat, RK <: ReproducingKernel_1}
Prepare and construct the spline.
Arguments
nodes: The function value nodes. This should be ann×n_1matrix, wherenis dimension of the sampled space andn_1is the number of function value nodes. It means that each column in the matrix defines one node.values: function values atnodesnodes.d_nodes: The function directional derivative nodes. This should be ann×n_2matrix, wherenis dimension of the sampled space andn_2is the number of function directional derivative nodes.es: Directions of the function directional derivatives. This should be ann×n_2matrix, wherenis dimension of the sampled space andn_2is the number of function directional derivative nodes. It means that each column in the matrix defines one direction of the function directional derivative.d_values: function directional derivative values atd_nodesnodes.kernel: reproducing kernel of Bessel potential space the normal spline is constructed in. It must be a struct object of the following type:RK_H1if the spline is constructing as a differentiable function,RK_H2if the spline is constructing as a twice differentiable function.
Return: the completely initialized NormalSpline object that can be passed to evaluate function.
interpolate(nodes::Vector{T}, values::Vector{T}, kernel::RK = RK_H0()) where {T <: AbstractFloat, RK <: ReproducingKernel_0}
Prepare and construct the 1D spline.
Arguments
nodes: The function value nodes. This should be ann×n_1matrix, wherenis dimension of the sampled space andn_1is the number of function value nodes. It means that each column in the matrix defines one node.values: function values atn_1interpolation nodes.kernel: reproducing kernel of Bessel potential space the normal spline is constructed in. It must be a struct object of the following type:RK_H0if the spline is constructing as a continuous function,RK_H1if the spline is constructing as a differentiable function,RK_H2if the spline is constructing as a twice differentiable function.
Return: the completely initialized NormalSpline object that can be passed to evaluate function.
interpolate(nodes::Vector{T}, values::Vector{T}, d_nodes::Vector{T}, d_values::Vector{T}, kernel::RK = RK_H1()) where {T <: AbstractFloat, RK <: ReproducingKernel_1}
Prepare and construct the 1D spline.
Arguments
nodes: The function value nodes.values: function values atnodesnodes.d_nodes: The function derivative nodes.d_values: function derivative values atd_nodesnodes.kernel: reproducing kernel of Bessel potential space the normal spline is constructed in. It must be a struct object of the following type:RK_H1if the spline is constructing as a differentiable function,RK_H2if the spline is constructing as a twice differentiable function.
Return: the completely initialized NormalSpline object that can be passed to evaluate function.
NormalSmoothingSplines.evaluate — Functionevaluate(spline::NormalSpline{T, RK}, points::Matrix{T}) where {T <: AbstractFloat, RK <: ReproducingKernel_0}
Evaluate the spline values at the locations defined in points.
Arguments
spline: theNormalSplineobject returned byinterpolateorconstruct` function.points: locations at which spline values are evaluating. This should be ann×mmatrix, wherenis dimension of the sampled space andmis the number of locations where spline values are evaluating. It means that each column in the matrix defines one location.
Return: Vector{T} of the spline values at the locations defined in points.
evaluate(spline::NormalSpline{T, RK}, points::Vector{T}) where {T <: AbstractFloat, RK <: ReproducingKernel_0}
Evaluate the spline values/value at the points locations.
Arguments
spline: theNormalSplineobject returned byinterpolateorconstructfunction.points: locations at which spline values are evaluating. This should be a vector of sizemwheremis the number of evaluating points.
Return: spline value at the point location.
NormalSmoothingSplines.evaluate_one — Functionevaluate_one(spline::NormalSpline{T, RK}, point::Vector{T}) where {T <: AbstractFloat, RK <: ReproducingKernel_0}
Evaluate the spline value at the point location.
Arguments
spline: theNormalSplineobject returned byinterpolateorconstructfunction.point: location at which spline value is evaluating. This should be a vector of sizen, wherenis dimension of the sampled space.
Return: the spline value at the location defined in point.
evaluate_one(spline::NormalSpline{T, RK}, point::T) where {T <: AbstractFloat, RK <: ReproducingKernel_0}
Evaluate the 1D spline value at the point location.
Arguments
spline: theNormalSplineobject returned byinterpolateorconstructfunction.point: location at which spline value is evaluating.
Return: spline value at the point location.
NormalSmoothingSplines.evaluate_gradient — Functionevaluate_gradient(spline::NormalSpline{T, RK}, point::Vector{T}) where {T <: AbstractFloat, RK <: ReproducingKernel_0}
Evaluate gradient of the spline at the location defined in point.
Arguments
spline: theNormalSplineobject returned byinterpolateorconstructfunction.point: location at which gradient value is evaluating. This should be a vector of sizen, wherenis dimension of the sampled space.
Note: Gradient of spline built with reproducing kernel RK_H0 does not exist at the spline nodes.
Return: Vector{T} - gradient of the spline at the location defined in point.
NormalSmoothingSplines.evaluate_derivative — Functionevaluate_derivative(spline::NormalSpline{T, RK}, point::T) where {T <: AbstractFloat, RK <: ReproducingKernel_0}
Evaluate the 1D spline derivative at the point location.
Arguments
spline: theNormalSplineobject returned byinterpolateorconstructfunction.point: location at which spline derivative is evaluating.
Note: Derivative of spline built with reproducing kernel RK_H0 does not exist at the spline nodes.
Return: the spline derivative value at the point location.
NormalSmoothingSplines.estimate_accuracy — Functionestimate_accuracy(spline::NormalSpline{T, RK}) where {T <: AbstractFloat, RK <: ReproducingKernel_0}
Assess accuracy of interpolation results by analyzing residuals.
Arguments
spline: theNormalSplineobject returned byconstructorinterpolatefunction.
Return: an estimation of the number of significant digits in the interpolation result.
NormalSmoothingSplines.estimate_cond — Functionestimate_cond(spline::NormalSpline{T, RK}) where {T <: AbstractFloat, RK <: ReproducingKernel_0}
Get an estimation of the Gram matrix condition number. It needs the spline object is prepared and requires O(N^2) operations. (C. Brás, W. Hager, J. Júdice, An investigation of feasible descent algorithms for estimating the condition number of a matrix. TOP Vol.20, No.3, 2012.)
Arguments
spline: theNormalSplineobject returned byprepare,constructorinterpolatefunction.
Return: an estimation of the Gram matrix condition number.
NormalSmoothingSplines.get_epsilon — Functionget_epsilon(spline::NormalSpline{T, RK}) where {T <: AbstractFloat, RK <: ReproducingKernel_0}
Get the 'scaling parameter' of Bessel Potential space the spline was built in.
Arguments
spline: theNormalSplineobject returned byprepare,constructorinterpolatefunction.
Return: ε - the 'scaling parameter'.
NormalSmoothingSplines.estimate_epsilon — Functionestimate_epsilon(nodes::Matrix{T}, kernel::RK = RK_H0()) where {T <: AbstractFloat, RK <: ReproducingKernel_0}
Get the estimation of the 'scaling parameter' of Bessel Potential space the spline being built in. It coincides with the result returned by get_epsilon function.
Arguments
nodes: The function value nodes. This should be ann×n_1matrix, wherenis dimension of the sampled space andn_1is the number of function value nodes. It means that each column in the matrix defines one node.kernel: reproducing kernel of Bessel potential space the normal spline will be constructed in. It must be a struct object of the following type:RK_H0if the spline is constructing as a continuous function,RK_H1if the spline is constructing as a differentiable function,RK_H2if the spline is constructing as a twice differentiable function.
Return: estimation of ε.
estimate_epsilon(nodes::Matrix{T}, d_nodes::Matrix{T}, kernel::RK = RK_H1()) where {T <: AbstractFloat, RK <: ReproducingKernel_1}
Get an the estimation of the 'scaling parameter' of Bessel Potential space the spline being built in. It coincides with the result returned by get_epsilon function.
Arguments
nodes: The function value nodes. This should be ann×n_1matrix, wherenis dimension of the sampled space andn_1is the number of function value nodes. It means that each column in the matrix defines one node.d_nodes: The function directional derivative nodes. This should be ann×n_2matrix, wherenis dimension of the sampled space andn_2is the number of function directional derivative nodes.kernel: reproducing kernel of Bessel potential space the normal spline will be constructed in. It must be a struct object of the following type:RK_H1if the spline is constructing as a differentiable function,RK_H2if the spline is constructing as a twice differentiable function.
Return: estimation of ε.
estimate_epsilon(nodes::Vector{T}, kernel::RK = RK_H0()) where {T <: AbstractFloat, RK <: ReproducingKernel_0}
Get an the estimation of the 'scaling parameter' of Bessel Potential space the spline being built in. It coincides with the result returned by get_epsilon function.
Arguments
nodes: The function value nodes.kernel: reproducing kernel of Bessel potential space the normal spline is constructed in. It must be a struct object of the following type:RK_H0if the spline is constructing as a continuous function,RK_H1if the spline is constructing as a differentiable function,RK_H2if the spline is constructing as a twice differentiable function.
Return: estimation of ε.
estimate_epsilon(nodes::Vector{T}, d_nodes::Vector{T}, kernel::RK = RK_H1()) where {T <: AbstractFloat, RK <: ReproducingKernel_1}
Get an the estimation of the 'scaling parameter' of Bessel Potential space the spline being built in. It coincides with the result returned by get_epsilon function.
Arguments
nodes: The function value nodes.d_nodes: The function derivative nodes.kernel: reproducing kernel of Bessel potential space the normal spline is constructed in. It must be a struct object of the following type:RK_H1if the spline is constructing as a differentiable function,RK_H2if the spline is constructing as a twice differentiable function.
Return: estimation of ε.
NormalSmoothingSplines.get_cond — Functionget_cond(nodes::Matrix{T}, kernel::RK = RK_H0()) where {T <: AbstractFloat, RK <: ReproducingKernel_0}
Get a value of the Gram matrix spectral condition number. It is obtained by means of the matrix SVD decomposition and requires $O(N^3)$ operations.
Arguments
nodes: The function value nodes. This should be ann×n_1matrix, wherenis dimension of the sampled space andn_1is the number of function value nodes. It means that each column in the matrix defines one node.kernel: reproducing kernel of Bessel potential space the normal spline is constructed in. It must be a struct object of the following type:RK_H0if the spline is constructing as a continuous function,RK_H1if the spline is constructing as a differentiable function,RK_H2if the spline is constructing as a twice differentiable function.
Return: a value of the Gram matrix spectral condition number.
get_cond(nodes::Matrix{T}, d_nodes::Matrix{T}, es::Matrix{T}, kernel::RK = RK_H1()) where {T <: AbstractFloat, RK <: ReproducingKernel_1}
Get a value of the Gram matrix spectral condition number. It is obtained by means of the matrix SVD decomposition and requires $O(N^3)$ operations.
Arguments
nodes: The function value nodes. This should be ann×n_1matrix, wherenis dimension of the sampled space andn_1is the number of function value nodes. It means that each column in the matrix defines one node.d_nodes: The function directional derivatives nodes. This should be ann×n_2matrix, wherenis dimension of the sampled space andn_2is the number of function directional derivative nodes.es: Directions of the function directional derivatives. This should be ann×n_2matrix, wherenis dimension of the sampled space andn_2is the number of function directional derivative nodes. It means that each column in the matrix defines one direction of the function directional derivative.kernel: reproducing kernel of Bessel potential space the normal spline is constructed in. It must be a struct object of the following type:RK_H1if the spline is constructing as a differentiable function,RK_H2if the spline is constructing as a twice differentiable function.
Return: a value of the Gram matrix spectral condition number.
Types
Bessel potential Space Reproducing Kernels
NormalSmoothingSplines.RK_H0 — Typestruct RK_H0{T <: AbstractFloat} <: ReproducingKernel_0
Defines a type of reproducing kernel of Bessel Potential space $H^{n/2 + 1/2}_ε (R^n)$ ('Basic Matérn kernel'):
Fields
ε::T: 'scaling parameter' from the Bessel Potential space definition, it may be omitted in the struct constructor otherwise it must be greater than zero
NormalSmoothingSplines.RK_H1 — Typestruct RK_H1{T <: AbstractFloat} <: ReproducingKernel_1
Defines a type of reproducing kernel of Bessel Potential space $H^{n/2 + 3/2}_ε (R^n)$ ('Linear Matérn kernel'):
Fields
ε::T: 'scaling parameter' from the Bessel Potential space definition, it may be omitted in the struct constructor otherwise it must be greater than zero
NormalSmoothingSplines.RK_H2 — Typestruct RK_H2{T <: AbstractFloat} <: ReproducingKernel_2
Defines a type of reproducing kernel of Bessel Potential space $H^{n/2 + 5/2}_ε (R^n)$ ('Quadratic Matérn kernel'):
Fields
ε::T: 'scaling parameter' from the Bessel Potential space definition, it may be omitted in the struct constructor otherwise it must be greater than zero
NormalSpline structure
NormalSmoothingSplines.NormalSpline — Typestruct NormalSpline{T, RK} <: AbstractSpline where {T <: AbstractFloat, RK <: ReproducingKernel_0}
Define a structure containing full information of a normal spline
Fields
_kernel: a reproducing kernel spline was built with_compression: factor of transforming the original node locations into unit hypercube_nodes: transformed function value nodes_values: function values at interpolation nodes_d_nodes: transformed function directional derivative nodes_es: normalized derivative directions_d_values: function directional derivative values_min_bound: minimal bounds of the original node locations area_gram: Gram matrix of the problem_chol: Cholesky factorization of the Gram matrix_mu: spline coefficients_cond: estimation of the Gram matrix condition number
Index
NormalSmoothingSplines.constructNormalSmoothingSplines.estimate_accuracyNormalSmoothingSplines.estimate_condNormalSmoothingSplines.estimate_epsilonNormalSmoothingSplines.evaluateNormalSmoothingSplines.evaluate_derivativeNormalSmoothingSplines.evaluate_gradientNormalSmoothingSplines.evaluate_oneNormalSmoothingSplines.get_condNormalSmoothingSplines.get_epsilonNormalSmoothingSplines.interpolateNormalSmoothingSplines.prepareNormalSmoothingSplines.NormalSplineNormalSmoothingSplines.RK_H0NormalSmoothingSplines.RK_H1NormalSmoothingSplines.RK_H2