Public API
NormalSplines.evaluate
NormalSplines.get_cond
NormalSplines.get_epsilon
NormalSplines.interpolate
NormalSplines.RK_H1
NormalSplines.RK_H2
NormalSplines.RK_H3
NormalSplines.RK_W1
NormalSplines.RK_W2
NormalSplines.RK_W3
Functions
NormalSplines.interpolate
— Functioninterpolate(knots::Vector{T}, values::Vector{T}, kernel :: RK, factorize = true) where {T <: AbstractFloat, RK <: ReproducingKernel_1}
Create an interpolating normal spline by values of function in knots.
Arguments
knots::Vector{T}
: locations of given function values.values::Vector{T}
: function values inknots
.kernel::RK
: reproducing kernel of Hilbert space the normal spline is constructed in. It must be a struct object of the following type:RK_W1
orRK_H1
if the spline is constructing as a continuous function,RK_W2
orRK_H2
if the spline is constructing as a differentiable function,RK_W3
orRK_H3
if the spline is constructing as a twice differentiable function.factorize::Bool
: defines if it is necessary to compute the Gram matrix factorization. Must be set totrue
if the spline is constructing first time. Can be set tofalse
if the spline is constructing with the same knots as it was done first time but with different function values in knots.
Returns : nothing.
interpolate(knots::Vector{T}, values::Vector{T}, d_knots::Vector{T}, d_values::Vector{T}, kernel :: RK, derivative = 1, factorize = true) where {T <: AbstractFloat, RK <: ReproducingKernel_2}
Create an interpolating normal spline by values of function and it's first or second derivatives in knots.
Arguments
knots::Vector{T}
: locations of given function values.values::Vector{T}
: function values inknots
.d_knots::Vector{T}
: locations of given function first derivative values (ifderivative
=1
) or second derivative values (ifderivative
=2
).d_values::Vector{T}
: values of function first derivative (if derivative=1
) or values of function second derivative (if derivative=2
) ind_knots
.kernel::RK
: reproducing kernel of Hilbert space the normal spline is constructed in. In a case whenderivative
=1
,kernel
must be a struct object of the following type:RK_W2
orRK_H2
if the spline is constructing as a differentiable function,RK_W3
orRK_H3
if the spline is constructing as a twice differentiable function. In a case whenderivative
=2
,kernel
must be a struct object of the following type:RK_W3
orRK_H3
- the spline is constructing as a twice differentiable function,derivative::Int
:1
-d_knots
andd_values
provide data of function first derivative,2
-d_knots
andd_values
provide data of function second derivative.factorize::Bool
: defines if it is necessary to compute the Gram matrix factorization. Must be set totrue
if the spline is constructing first time. Can be set tofalse
if the spline is constructing with the same knots as it was done first time but with different function values in knots.
Returns : nothing.
interpolate(knots::Vector{T}, values::Vector{T}, d1_knots::Vector{T}, d1_values::Vector{T}, d2_knots::Vector{T}, d2_values::Vector{T}, kernel :: RK, factorize = true) where {T <: AbstractFloat, RK <: ReproducingKernel_3}
Create an interpolating normal spline by values of function and it's first and second derivatives in knots.
Arguments
knots::Vector{T}
: locations of given function values.values::Vector{T}
: function values inknots
.d1_knots::Vector{T}
: locations of given function derivative values.d1_values::Vector{T}
: values of function derivative ind1_knots
.d2_knots::Vector{T}
: locations of given function second derivative values.d2_values::Vector{T}
: values of function second derivative ind2_knots
.kernel::RK
: reproducing kernel of Hilbert space the normal spline is constructed in. It must be a struct object of the following type:RK_W3
orRK_H3
- the spline is constructing as a twice differentiable function.factorize::Bool
: defines if it is necessary to compute the Gram matrix factorization. Must be set totrue
if the spline is constructing first time. Can be set tofalse
if the spline is constructing with the same knots as it was done first time but with different function values in knots.
Returns : nothing.
NormalSplines.evaluate
— Functionevaluate(points::Vector{T}, derivative::Int = 0) where T <: AbstractFloat
Evaluate normal spline and it's derivatives.
Arguments
points::Vector{T}
: locations at which spline or its derivative are evaluating.derivative::Int
:0
- spline evaluation.1
- first derivative of spline evaluation.2
- second derivative of spline evaluation.
Returns: Vector{T}
of spline values or its derivatives at the locations defined in points
.
evaluate(point::T, derivative::Int = 0) where T <: AbstractFloat
Evaluate normal spline and it's derivatives.
Arguments
point::T
: location at which spline or its derivative are evaluating.derivative::Int
:0
- spline evaluation.1
- first derivative of spline evaluation.2
- second derivative of spline evaluation.
Returns: spline value or its derivative at the point
location.
NormalSplines.get_cond
— Functionget_cond()
Estimate condition number of Gram matrix.
NormalSplines.get_epsilon
— Functionget_epsilon()
Get value of parameter ε
of struct object RK_H1
, RK_H2
or RK_H3
.
Types
Sobolev space Reproducing Kernels
NormalSplines.RK_W1
— Typestruct RK_W1 <: ReproducingKernel_1
Define a type of reproducing kernel of Sobolev space $W^2_1 [0,1]$:
NormalSplines.RK_W2
— Typestruct RK_W2 <: ReproducingKernel_2
Define a type of reproducing kernel of Sobolev space $W^2_2 [0,1]$:
NormalSplines.RK_W3
— Typestruct RK_W3 <: ReproducingKernel_3
Define a type of reproducing kernel of Sobolev space $W^2_3 [0,1]$:
Bessel Potential space Reproducing Kernels
NormalSplines.RK_H1
— Typestruct RK_H1{T <: AbstractFloat} <: ReproducingKernel_1
Define a type of reproducing kernel of Bessel Potential space $H^1_ε (R)$:
Arguments
- no argument: means that parameter
ε
will be estimated in course of interpolation procedure execution. ε::T
: 'scaling parameter' from the Bessel Potential space definition, it must be greater than zero.
NormalSplines.RK_H2
— Typestruct RK_H2{T <: AbstractFloat} <: ReproducingKernel_2
Define a type of reproducing kernel of Bessel Potential space $H^2_ε (R)$:
Arguments
- no argument: means that parameter
ε
will be estimated in course of interpolation procedure execution. ε::T
: 'scaling parameter' from the Bessel Potential space definition, it must be greater than zero.
NormalSplines.RK_H3
— Typestruct RK_H3{T <: AbstractFloat} <: ReproducingKernel_3
Define a type of reproducing kernel of Bessel Potential space $H^3_ε (R)$:
Arguments
- no argument: means that parameter
ε
will be estimated in course of interpolation procedure execution. ε::T
: 'scaling parameter' from the Bessel Potential space definition, it must be greater than zero.