Public API

API Summary

FunctionDescription
prepare_approximationPrepare the smoothing spline by constructing and factoring a Gram matrix of the problem.
construct_approximationConstruct the smoothing spline by calculating its coefficients.
approximatePrepare and construct the smoothing spline.
prepare_interpolationPrepare the interpolating spline by constructing and factoring a Gram matrix of the problem.
construct_interpolationConstruct the interpolating spline by calculating its coefficients.
interpolatePrepare and construct the interpolating spline.
evaluateEvaluate the spline values.
evaluate_atEvaluate the spline value at the required location.
evaluate_gradientEvaluate gradient of the spline at the required location.
evaluate_derivativeEvaluate the 1-D spline derivative at the required location.
estimate_accuracyEstimate accuracy of the function interpolation result.
estimate_condEstimate the Gram matrix 1-norm condition number.
estimate_epsilonEstimate the 'scaling parameter' of Bessel potential space the spline being built in.
get_epsilonGet the 'scaling parameter' of Bessel potential space the normal spline was built in.
get_condGet the Gram matrix spectral condition number.

Functions

NormalSmoothingSplines.prepare_approximationFunction

prepare_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 an n×n_1_b matrix, where n is dimension of the sampled space and n_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.

source

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 an n×n_1 matrix, where n is dimension of the sampled space and n_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 an n×n_1_b matrix, where n is dimension of the sampled space and n_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.

source

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 an n_1_b vector where n_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.

source

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 an n_1 vector where n_1 is the number of function value interpolation nodes.
  • nodes_b: function value approximation nodes. This should be an n_1_b vector where n_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.

source
NormalSmoothingSplines.construct_approximationFunction

construct_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 initialized NormalSpline object returned by prepare function.
  • values_lb: function lower bound values at approximation nodes
  • values_ub: function upper bound values at approximation nodes
  • 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 and the number of QP algorithm iterations done.

source

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 initialized NormalSpline object returned by prepare function.
  • values: function values at interpolation nodes.
  • values_lb: function lower bound values at approximation nodes
  • values_ub: function upper bound values at approximation nodes
  • 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 and the number of QP algorithm iterations done.

source
NormalSmoothingSplines.approximateFunction

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_b: function value approximation nodes. This should be an n×n_1_b matrix, where n is dimension of the sampled space and n_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 nodes
  • values_ub: function upper bound values at 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.
  • 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.

source

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 an n×n_1 matrix, where n is dimension of the sampled space and n_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 an n×n_1_b matrix, where n is dimension of the sampled space and n_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 nodes
  • values_ub: function upper bound values at 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.
  • 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.

source

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 an n_1_b vector where n_1_b is the number of function value approximation nodes.
  • values_lb: function lower bound values at approximation nodes
  • values_ub: function upper bound values at 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.
  • 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.

source

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 an n_1 vector where n_1 is the number of function value interpolation nodes.
  • values: function values at interpolation nodes.
  • nodes_b: function value approximation nodes. This should be an n_1_b vector where n_1_b is the number of function value approximation nodes.
  • values_lb: function lower bound values at approximation nodes
  • values_ub: function upper bound values at 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.
  • 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.

source
NormalSmoothingSplines.prepare_interpolationFunction

prepare_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 an n×n_1 matrix, where n is dimension of the sampled space and n_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.

source

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 an n×n_1 matrix, where n is dimension of the sampled space and n_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 an n×n_2 matrix, where n is dimension of the sampled space and n_2 is the number of function directional derivative nodes.
  • es: Directions of the function directional derivatives. This should be an n×n_2 matrix, where n is dimension of the sampled space and n_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.

source

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 an n_1 vector where n_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.

source

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 an n_1 vector where n_1 is the number of function value nodes.
  • d_nodes: The function derivatives nodes. This should be an n_2 vector where n_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.

source
NormalSmoothingSplines.construct_interpolationFunction

construct_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 prepared NormalSpline object.
  • values: function values at interpolation nodes.

Return: constructed NormalSpline object.

source

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 prepared NormalSpline object.
  • values: function values at interpolation nodes.
  • d_values: function directional derivative values at function derivative nodes.

Return: constructed NormalSpline object.

source
NormalSmoothingSplines.interpolateFunction

interpolate(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 an n×n_1 matrix, where n is dimension of the sampled space and n_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.

source

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 an n×n_1 matrix, where n is dimension of the sampled space and n_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 an n×n_2 matrix, where n is dimension of the sampled space and n_2 is the number of function directional derivative nodes.
  • es: Directions of the function directional derivatives. This should be an n×n_2 matrix, where n is dimension of the sampled space and n_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.

source

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 an n_1 vector where n_1 is the number of function value nodes.
  • values: function values at n_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.

source

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 an n_1 vector where n_1 is the number of function value nodes.
  • values: function values at nodes nodes.
  • d_nodes: The function derivatives nodes. This should be an n_2 vector where n_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.

source
NormalSmoothingSplines.evaluateFunction

evaluate(spline::NormalSpline{T, RK}, points::Matrix{T}) where {T <: AbstractFloat, RK <: ReproducingKernel_0}

Evaluate the spline values at points locations.

Arguments

  • spline: constructedNormalSpline` object.
  • points: locations at which spline values are evaluating. This should be an n×m matrix, where n is dimension of the sampled space and m 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.

source

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: constructed NormalSpline object.
  • points: locations at which spline values are evaluating. This should be a vector of size m where m is the number of evaluating points.

Return: spline value at the point location.

source
NormalSmoothingSplines.evaluate_atFunction

evaluate_at(spline::NormalSpline{T, RK}, point::Vector{T}) where {T <: AbstractFloat, RK <: ReproducingKernel_0}

Evaluate the spline value at the point location.

Arguments

  • spline: constructed NormalSpline object.
  • point: location at which spline value is evaluating. This should be a vector of size n, where n is dimension of the sampled space.

Return: spline value at the location defined in point.

source

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: constructed NormalSpline object.
  • point: location at which spline value is evaluating.

Return: spline value at the point location.

source
NormalSmoothingSplines.evaluate_gradientFunction

evaluate_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: the NormalSpline object returned by interpolate or construct function.
  • point: location at which gradient value is evaluating. This should be a vector of size n, where n 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.

source
NormalSmoothingSplines.evaluate_derivativeFunction

evaluate_derivative(spline::NormalSpline{T, RK}, point::T) where {T <: AbstractFloat, RK <: ReproducingKernel_0}

Evaluate the 1D spline derivative at the point location.

Arguments

  • spline: the NormalSpline object returned by interpolate or construct 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.

source
NormalSmoothingSplines.estimate_accuracyFunction

estimate_accuracy(spline::NormalSpline{T, RK}) where {T <: AbstractFloat, RK <: ReproducingKernel_0}

Assess accuracy of interpolation results by analyzing residuals.

Arguments

  • spline: constructed NormalSpline object.

Return: estimation of the number of significant digits in the interpolation result.

source
NormalSmoothingSplines.estimate_condFunction

estimate_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: prepared NormalSpline object.

Return: an estimation of the Gram matrix condition number.

source
NormalSmoothingSplines.get_epsilonFunction

get_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: prepared NormalSpline object.

Return: ε - the 'scaling parameter'.

source
NormalSmoothingSplines.estimate_epsilonFunction

estimate_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 an n×n_1 matrix, where n is dimension of the sampled space and n_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 ε.

source

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 an n×n_1 matrix, where n is dimension of the sampled space and n_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 an n×n_2 matrix, where n is dimension of the sampled space and n_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 ε.

source

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 ε.

source

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 ε.

source
NormalSmoothingSplines.get_condFunction

get_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 an n×n_1 matrix, where n is dimension of the sampled space and n_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.

source

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 an n×n_1 matrix, where n is dimension of the sampled space and n_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 an n×n_2 matrix, where n is dimension of the sampled space and n_2 is the number of function directional derivative nodes.
  • es: Directions of the function directional derivatives. This should be an n×n_2 matrix, where n is dimension of the sampled space and n_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.

source

Types

Bessel potential Space Reproducing Kernels

NormalSmoothingSplines.RK_H0Type

struct 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'):

\[V(\eta , \xi, \varepsilon) = \exp (-\varepsilon |\xi - \eta|) \, .\]

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
source
NormalSmoothingSplines.RK_H1Type

struct 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'):

\[V(\eta , \xi, \varepsilon) = \exp (-\varepsilon |\xi - \eta|) (1 + \varepsilon |\xi - \eta|) \, .\]

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
source
NormalSmoothingSplines.RK_H2Type

struct 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'):

\[V(\eta , \xi, \varepsilon) = \exp (-\varepsilon |\xi - \eta|) (3 + 3\varepsilon |\xi - \eta| + \varepsilon ^2 |\xi - \eta| ^2 ) \, .\]

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
source

NormalSpline structure

NormalSmoothingSplines.NormalSplineType

struct 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.
source

Index