Public API
API Summary
Function | Description |
---|---|
prepare_approximation | Prepare the smoothing spline by constructing and factoring a Gram matrix of the problem. |
construct_approximation | Construct the smoothing spline by calculating its coefficients. |
approximate | Prepare and construct the smoothing spline. |
prepare_interpolation | Prepare the interpolating spline by constructing and factoring a Gram matrix of the problem. |
construct_interpolation | Construct the interpolating spline by calculating its coefficients. |
interpolate | Prepare and construct the interpolating spline. |
evaluate | Evaluate the spline values. |
evaluate_at | Evaluate the spline value at the required location. |
evaluate_gradient | Evaluate gradient of the spline at the required location. |
evaluate_derivative | Evaluate the 1-D 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_approximation
— Functionprepare_approximation(nodes_b::Matrix{T}, kernel::RK = RK_H0()) where {T <: AbstractFloat, RK <: ReproducingKernel_0}
Prepare the approximating normal spline by constructing and factoring a Gram matrix of the problem. Initialize the NormalSpline
object.
Arguments
nodes_b
: function value approximation nodes. This should be ann×n_1_b
matrix, wheren
is dimension of the sampled space andn_1_b
is the number of function value approximation 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_H0
if the spline is constructing as a continuous function,RK_H1
if the spline is constructing as a differentiable function,RK_H2
if the spline is constructing as a twice differentiable function.
Return: prepared NormalSpline
object.
prepare_approximation(nodes::Matrix{T}, nodes_b::Matrix{T}, kernel::RK = RK_H0()) where {T <: AbstractFloat, RK <: ReproducingKernel_0}
Prepare the approximating normal spline by constructing and factoring a Gram matrix of the problem. Initialize the NormalSpline
object.
Arguments
nodes
: function value interpolation nodes. This should be ann×n_1
matrix, wheren
is dimension of the sampled space andn_1
is the number of function value interpolation nodes. It means that each column in the matrix defines one node.nodes_b
: function value approximation nodes. This should be ann×n_1_b
matrix, wheren
is dimension of the sampled space andn_1_b
is the number of function value approximation 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_H0
if the spline is constructing as a continuous function,RK_H1
if the spline is constructing as a differentiable function,RK_H2
if the spline is constructing as a twice differentiable function.
Return: prepared NormalSpline
object.
prepare_approximation(nodes_b::Vector{T}, kernel::RK = RK_H0()) where {T <: AbstractFloat, RK <: ReproducingKernel_0}
Prepare the 1D approximating normal spline by constructing and factoring a Gram matrix of the problem. Initialize the NormalSpline
object.
Arguments
nodes_b
: function value approximation nodes. This should be ann_1_b
vector wheren_1_b
is the number of function value approximation 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_H0
if the spline is constructing as a continuous function,RK_H1
if the spline is constructing as a differentiable function,RK_H2
if the spline is constructing as a twice differentiable function.
Return: prepared NormalSpline
object.
prepare_approximation(nodes::Vector{T}, nodes_b::Vector{T}, kernel::RK = RK_H0()) where {T <: AbstractFloat, RK <: ReproducingKernel_0}
Prepare the 1D approximating normal spline by constructing and factoring a Gram matrix of the problem. Initialize the NormalSpline
object.
Arguments
nodes
: function value interpolation nodes. This should be ann_1
vector wheren_1
is the number of function value interpolation nodes.nodes_b
: function value approximation nodes. This should be ann_1_b
vector wheren_1_b
is the number of function value approximation 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_H0
if the spline is constructing as a continuous function,RK_H1
if the spline is constructing as a differentiable function,RK_H2
if the spline is constructing as a twice differentiable function.
Return: prepared NormalSpline
object.
NormalSmoothingSplines.construct_approximation
— Functionconstruct_approximation(spline::NormalSpline{T, RK}, values_lb::Vector{T}, values_ub::Vector{T}, maxiter::Int, ftol::T) where {T <: AbstractFloat, RK <: ReproducingKernel_0}
construct the approximating normal spline by calculating its coefficients and completely initializing the NormalSpline
object.
Arguments
spline
: the partly initializedNormalSpline
object returned byprepare
function.values_lb
: function lower bound values at approximation nodesvalues_ub
: function upper bound values at approximation nodesmaxiter
: Maximum allowed number of iterations.ftol
: convergence tolerance. The iteration stops when relative spline norm change is smaller than ftol.
Return: constructed NormalSpline
object and the number of QP algorithm iterations done.
construct_approximation(spline::NormalSpline{T, RK}, values::Vector{T}, values_lb::Vector{T}, values_ub::Vector{T}, maxiter::Int, ftol::T) where {T <: AbstractFloat, RK <: ReproducingKernel_0}
construct the approximating normal spline by calculating its coefficients and completely initializing the NormalSpline
object.
Arguments
spline
: the partly initializedNormalSpline
object returned byprepare
function.values
: function values at interpolation nodes.values_lb
: function lower bound values at approximation nodesvalues_ub
: function upper bound values at approximation nodesmaxiter
: Maximum allowed number of iterations.ftol
: convergence tolerance. The iteration stops when relative spline norm change is smaller than ftol.
Return: constructed NormalSpline
object and the number of QP algorithm iterations done.
NormalSmoothingSplines.approximate
— Functionapproximate(nodes::Matrix{T}, values::Vector{T}, nodes_b::Matrix{T}, values_lb::Vector{T}, values_ub::Vector{T}, kernel::RK, maxiter::Int, ftol::T) where {T <: AbstractFloat, RK <: ReproducingKernel_0}
Prepare and construct the approximating normal spline.
Arguments
nodes_b
: function value approximation nodes. This should be ann×n_1_b
matrix, wheren
is dimension of the sampled space andn_1_b
is the number of function value approximation nodes. It means that each column in the matrix defines one node.values_lb
: function lower bound values at approximation nodesvalues_ub
: function upper bound values at approximation nodeskernel
: reproducing kernel of Bessel potential space the normal spline is constructed in. It must be a struct object of the following type:RK_H0
if the spline is constructing as a continuous function,RK_H1
if the spline is constructing as a differentiable function,RK_H2
if the spline is constructing as a twice differentiable function.maxiter
: Maximum allowed number of iterations.ftol
: convergence tolerance. The iteration stops when relative spline norm change is smaller than ftol.
Return: constructed NormalSpline
object.
approximate(nodes::Matrix{T}, values::Vector{T}, nodes_b::Matrix{T}, values_lb::Vector{T}, values_ub::Vector{T}, kernel::RK, maxiter::Int, ftol::T) where {T <: AbstractFloat, RK <: ReproducingKernel_0}
Prepare and construct the approximating normal spline.
Arguments
nodes
: function value interpolation nodes. This should be ann×n_1
matrix, wheren
is dimension of the sampled space andn_1
is the number of function value nodes. It means that each column in the matrix defines one node.nodes_b
: function value approximation nodes. This should be ann×n_1_b
matrix, wheren
is dimension of the sampled space andn_1_b
is the number of function value approximation nodes. It means that each column in the matrix defines one node.values
: function values at interpolation nodes.values_lb
: function lower bound values at approximation nodesvalues_ub
: function upper bound values at approximation nodeskernel
: reproducing kernel of Bessel potential space the normal spline is constructed in. It must be a struct object of the following type:RK_H0
if the spline is constructing as a continuous function,RK_H1
if the spline is constructing as a differentiable function,RK_H2
if the spline is constructing as a twice differentiable function.maxiter
: Maximum allowed number of iterations.ftol
: convergence tolerance. The iteration stops when relative spline norm change is smaller than ftol.
Return: constructed NormalSpline
object.
approximate(values::Vector{T}, nodes_b::Vector{T}, values_lb::Vector{T}, values_ub::Vector{T}, kernel::RK, maxiter::Int, ftol::T) where {T <: AbstractFloat, RK <: ReproducingKernel_0}
Prepare and construct the 1D approximating normal spline.
Arguments
nodes_b
: function value approximation nodes. This should be ann_1_b
vector wheren_1_b
is the number of function value approximation nodes.values_lb
: function lower bound values at approximation nodesvalues_ub
: function upper bound values at approximation nodeskernel
: reproducing kernel of Bessel potential space the normal spline is constructed in. It must be a struct object of the following type:RK_H0
if the spline is constructing as a continuous function,RK_H1
if the spline is constructing as a differentiable function,RK_H2
if the spline is constructing as a twice differentiable function.maxiter
: Maximum allowed number of iterations.ftol
: convergence tolerance. The iteration stops when relative spline norm change is smaller than ftol.
Return: constructed NormalSpline
object.
approximate(nodes::Vector{T}, values::Vector{T}, nodes_b::Vector{T}, values_lb::Vector{T}, values_ub::Vector{T}, kernel::RK, maxiter::Int, ftol::T) where {T <: AbstractFloat, RK <: ReproducingKernel_0}
Prepare and construct the 1D approximating normal spline.
Arguments
nodes
: function value interpolation nodes. This should be ann_1
vector wheren_1
is the number of function value interpolation nodes.values
: function values at interpolation nodes.nodes_b
: function value approximation nodes. This should be ann_1_b
vector wheren_1_b
is the number of function value approximation nodes.values_lb
: function lower bound values at approximation nodesvalues_ub
: function upper bound values at approximation nodeskernel
: reproducing kernel of Bessel potential space the normal spline is constructed in. It must be a struct object of the following type:RK_H0
if the spline is constructing as a continuous function,RK_H1
if the spline is constructing as a differentiable function,RK_H2
if the spline is constructing as a twice differentiable function.maxiter
: Maximum allowed number of iterations.ftol
: convergence tolerance. The iteration stops when relative spline norm change is smaller than ftol.
Return: constructed NormalSpline
object.
NormalSmoothingSplines.prepare_interpolation
— Functionprepare_interpolation(nodes::Matrix{T}, kernel::RK = RK_H0()) where {T <: AbstractFloat, RK <: ReproducingKernel_0}
Prepare the interpolating normal spline by constructing and factoring Gram matrix of the problem. Initialize the NormalSpline
object.
Arguments
nodes
: function value interpolation nodes. This should be ann×n_1
matrix, wheren
is dimension of the sampled space andn_1
is 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_H0
if the spline is constructing as a continuous function,RK_H1
if the spline is constructing as a differentiable function,RK_H2
if the spline is constructing as a twice differentiable function.
Return: prepared NormalSpline
object.
prepare_interpolation(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
: function value interpolation nodes. This should be ann×n_1
matrix, wheren
is dimension of the sampled space andn_1
is the number of function value nodes. It means that each column in the matrix defines one node.d_nodes
: function directional derivatives nodes. This should be ann×n_2
matrix, wheren
is dimension of the sampled space andn_2
is the number of function directional derivative nodes.es
: Directions of the function directional derivatives. This should be ann×n_2
matrix, wheren
is dimension of the sampled space andn_2
is 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_H1
if the spline is constructing as a differentiable function,RK_H2
if the spline is constructing as a twice differentiable function.
Return: prepared NormalSpline
object.
prepare_interpolation(nodes::Vector{T}, kernel::RK = RK_H0()) where {T <: AbstractFloat, RK <: ReproducingKernel_0}
Prepare the 1D interpolating normal spline by constructing and factoring a Gram matrix of the problem. Initialize the NormalSpline
object.
Arguments
nodes
: function value interpolation nodes. This should be ann_1
vector wheren_1
is the number of 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_H0
if the spline is constructing as a continuous function,RK_H1
if the spline is constructing as a differentiable function,RK_H2
if the spline is constructing as a twice differentiable function.
Return: prepared NormalSpline
object.
prepare_interpolation(nodes::Vector{T}, d_nodes::Vector{T}, kernel::RK = RK_H1()) where {T <: AbstractFloat, RK <: ReproducingKernel_1}
Prepare the 1D interpolating normal spline by constructing and factoring a Gram matrix of the problem. Initialize the NormalSpline
object.
Arguments
nodes
: function value interpolation nodes. This should be ann_1
vector wheren_1
is the number of function value nodes.d_nodes
: The function derivatives nodes. This should be ann_2
vector wheren_2
is the number of 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_H1
if the spline is constructing as a differentiable function,RK_H2
if the spline is constructing as a twice differentiable function.
Return: prepared NormalSpline
object.
NormalSmoothingSplines.construct_interpolation
— Functionconstruct_interpolation(spline::NormalSpline{T, RK}, values::Vector{T}) where {T <: AbstractFloat, RK <: ReproducingKernel_0}
Construct the interpolating normal spline by calculating its coefficients and completely initializing the NormalSpline
object.
Arguments
spline
: the preparedNormalSpline
object.values
: function values at interpolation nodes.
Return: constructed NormalSpline
object.
construct_interpolation(spline::NormalSpline{T, RK}, values::Vector{T}, d_values::Vector{T}) where {T <: AbstractFloat, RK <: ReproducingKernel_1}
Construct the interpolating normal spline by calculating its coefficients and completely initializing the NormalSpline
object.
Arguments
spline
: the preparedNormalSpline
object.values
: function values at interpolation nodes.d_values
: function directional derivative values at function derivative nodes.
Return: constructed NormalSpline
object.
NormalSmoothingSplines.interpolate
— Functioninterpolate(nodes::Matrix{T}, values::Vector{T}, kernel::RK = RK_H0()) where {T <: AbstractFloat, RK <: ReproducingKernel_0}
Prepare and construct the interpolating normal spline.
Arguments
nodes
: function value interpolation nodes. This should be ann×n_1
matrix, wheren
is dimension of the sampled space andn_1
is the number of function value nodes. It means that each column in the matrix defines one node.values
: function values at interpolation 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_H0
if the spline is constructing as a continuous function,RK_H1
if the spline is constructing as a differentiable function,RK_H2
if the spline is constructing as a twice differentiable function.
Return: constructed NormalSpline
object.
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 interpolating normal spline.
Arguments
nodes
: function value interpolation nodes. This should be ann×n_1
matrix, wheren
is dimension of the sampled space andn_1
is the number of function value nodes. It means that each column in the matrix defines one node.values
: function values at interpolation nodes.d_nodes
: function directional derivative nodes. This should be ann×n_2
matrix, wheren
is dimension of the sampled space andn_2
is the number of function directional derivative nodes.es
: Directions of the function directional derivatives. This should be ann×n_2
matrix, wheren
is dimension of the sampled space andn_2
is 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 at 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_H1
if the spline is constructing as a differentiable function,RK_H2
if the spline is constructing as a twice differentiable function.
Return: constructed NormalSpline
object.
interpolate(nodes::Vector{T}, values::Vector{T}, kernel::RK = RK_H0()) where {T <: AbstractFloat, RK <: ReproducingKernel_0}
Prepare and construct the 1D interpolating normal spline.
Arguments
nodes
: function value interpolation nodes. This should be ann_1
vector wheren_1
is the number of function value nodes.values
: function values atn_1
interpolation 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_H0
if the spline is constructing as a continuous function,RK_H1
if the spline is constructing as a differentiable function,RK_H2
if the spline is constructing as a twice differentiable function.
Return: constructed NormalSpline
object.
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_interpolation and construct the 1D interpolating normal spline.
Arguments
nodes
: function value interpolation nodes. This should be ann_1
vector wheren_1
is the number of function value nodes.values
: function values atnodes
nodes.d_nodes
: The function derivatives nodes. This should be ann_2
vector wheren_2
is the number of function derivatives nodes.d_values
: function derivative values at 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_H1
if the spline is constructing as a differentiable function,RK_H2
if the spline is constructing as a twice differentiable function.
Return: constructed NormalSpline
object.
NormalSmoothingSplines.evaluate
— Functionevaluate(spline::NormalSpline{T, RK}, points::Matrix{T}) where {T <: AbstractFloat, RK <: ReproducingKernel_0}
Evaluate the spline values at points
locations.
Arguments
spline: constructed
NormalSpline` object.points
: locations at which spline values are evaluating. This should be ann×m
matrix, wheren
is dimension of the sampled space andm
is 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 1D spline values at the points
locations.
Arguments
spline
: constructedNormalSpline
object.points
: locations at which spline values are evaluating. This should be a vector of sizem
wherem
is the number of evaluating points.
Return: spline value at the point
location.
NormalSmoothingSplines.evaluate_at
— Functionevaluate_at(spline::NormalSpline{T, RK}, point::Vector{T}) where {T <: AbstractFloat, RK <: ReproducingKernel_0}
Evaluate the spline value at the point
location.
Arguments
spline
: constructedNormalSpline
object.point
: location at which spline value is evaluating. This should be a vector of sizen
, wheren
is dimension of the sampled space.
Return: spline value at the location defined in point
.
evaluate_at(spline::NormalSpline{T, RK}, point::T) where {T <: AbstractFloat, RK <: ReproducingKernel_0}
Evaluate the 1D spline value at the point
location.
Arguments
spline
: constructedNormalSpline
object.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
: theNormalSpline
object returned byinterpolate
orconstruct
function.point
: location at which gradient value is evaluating. This should be a vector of sizen
, wheren
is 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 point
location.
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
: theNormalSpline
object returned byinterpolate
orconstruct
function.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: 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
: constructedNormalSpline
object.
Return: 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
: preparedNormalSpline
object.
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
: preparedNormalSpline
object.
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'2 of Bessel Potential space the spline being built in. It coincides with the result returned by get_epsilon
function.
Arguments
nodes
: function value interpolation nodes. This should be ann×n_1
matrix, wheren
is dimension of the sampled space andn_1
is 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_H0
if the spline is constructing as a continuous function,RK_H1
if the spline is constructing as a differentiable function,RK_H2
if 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
: function value interpolation nodes. This should be ann×n_1
matrix, wheren
is dimension of the sampled space andn_1
is the number of function value nodes. It means that each column in the matrix defines one node.d_nodes
: function directional derivative nodes. This should be ann×n_2
matrix, wheren
is dimension of the sampled space andn_2
is 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_H1
if the spline is constructing as a differentiable function,RK_H2
if 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
: function value interpolation 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_H0
if the spline is constructing as a continuous function,RK_H1
if the spline is constructing as a differentiable function,RK_H2
if 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
: function value interpolation nodes.d_nodes
: 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_H1
if the spline is constructing as a differentiable function,RK_H2
if 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.
Arguments
nodes
: function value interpolation nodes. This should be ann×n_1
matrix, wheren
is dimension of the sampled space andn_1
is 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_H0
if the spline is constructing as a continuous function,RK_H1
if the spline is constructing as a differentiable function,RK_H2
if 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.
Arguments
nodes
: function value interpolation nodes. This should be ann×n_1
matrix, wheren
is dimension of the sampled space andn_1
is the number of function value nodes. It means that each column in the matrix defines one node.d_nodes
: function directional derivatives nodes. This should be ann×n_2
matrix, wheren
is dimension of the sampled space andn_2
is the number of function directional derivative nodes.es
: Directions of the function directional derivatives. This should be ann×n_2
matrix, wheren
is dimension of the sampled space andn_2
is 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_H1
if the spline is constructing as a differentiable function,RK_H2
if the spline is constructing as a twice differentiable function.
Return: 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 interpolation nodes_nodes_b
: transformed function value approximation nodes_values
: function values at interpolation nodes_values_lb
: function lower bound values at approximation nodes_values_ub
: function upper bound values at approximation nodes_d_nodes
: transformed function directional derivative interpolation nodes_d_nodes_b
: transformed function directional derivative approximation nodes_es
: normalized derivative directions at interpolation nodes_es_b
: normalized derivative directions at approximation nodes_d_values
: function directional derivative values at interpolation nodes_d_values_lb
: function lower bound directional derivative values at approximation nodes_d_values_ub
: function upper bound directional derivative values at approximation nodes_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_active
: active inequality constraint numbers at solution_cond
: estimation of the Gram matrix condition number_ier
: An integer flag. If it is equal to 0, the optimal solution was found. If it is equal to 1, the approximate solution was found. QP algorithm iterations were stopped because of small spline norm change. If it is equal to 2, the approximate solution was found. QP algorithm iterations were stopped because of maximum allowed number of iterations was reached. If it is equal to -1, the solution was not calculated.
Index
NormalSmoothingSplines.approximate
NormalSmoothingSplines.construct_approximation
NormalSmoothingSplines.construct_interpolation
NormalSmoothingSplines.estimate_accuracy
NormalSmoothingSplines.estimate_cond
NormalSmoothingSplines.estimate_epsilon
NormalSmoothingSplines.evaluate
NormalSmoothingSplines.evaluate_at
NormalSmoothingSplines.evaluate_derivative
NormalSmoothingSplines.evaluate_gradient
NormalSmoothingSplines.get_cond
NormalSmoothingSplines.get_epsilon
NormalSmoothingSplines.interpolate
NormalSmoothingSplines.prepare_approximation
NormalSmoothingSplines.prepare_interpolation
NormalSmoothingSplines.NormalSpline
NormalSmoothingSplines.RK_H0
NormalSmoothingSplines.RK_H1
NormalSmoothingSplines.RK_H2