Public API

Functions

NormalSplines.interpolateFunction

interpolate(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 in 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_W1 or RK_H1 if the spline is constructing as a continuous function, RK_W2 or RK_H2 if the spline is constructing as a differentiable function, RK_W3 or RK_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 to true if the spline is constructing first time. Can be set to false if the spline is constructing with the same knots as it was done first time but with different function values in knots.

Returns : nothing.

source

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 in knots.
  • d_knots::Vector{T}: locations of given function first derivative values (if derivative=1) or second derivative values (if derivative=2).
  • d_values::Vector{T}: values of function first derivative (if derivative=1) or values of function second derivative (if derivative=2) in d_knots.
  • kernel::RK: reproducing kernel of Hilbert space the normal spline is constructed in. In a case when derivative=1, kernel must be a struct object of the following type: RK_W2 or RK_H2 if the spline is constructing as a differentiable function, RK_W3 or RK_H3 if the spline is constructing as a twice differentiable function. In a case when derivative=2, kernel must be a struct object of the following type: RK_W3 or RK_H3 - the spline is constructing as a twice differentiable function,
  • derivative::Int: 1 - d_knots and d_values provide data of function first derivative, 2 - d_knots and d_values provide data of function second derivative.
  • factorize::Bool: defines if it is necessary to compute the Gram matrix factorization. Must be set to true if the spline is constructing first time. Can be set to false if the spline is constructing with the same knots as it was done first time but with different function values in knots.

Returns : nothing.

source

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 in knots.
  • d1_knots::Vector{T}: locations of given function derivative values.
  • d1_values::Vector{T}: values of function derivative in d1_knots.
  • d2_knots::Vector{T}: locations of given function second derivative values.
  • d2_values::Vector{T}: values of function second derivative in d2_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 or RK_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 to true if the spline is constructing first time. Can be set to false if the spline is constructing with the same knots as it was done first time but with different function values in knots.

Returns : nothing.

source
NormalSplines.evaluateFunction

evaluate(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.

source

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.

source

Types

Sobolev space Reproducing Kernels

NormalSplines.RK_W1Type

struct RK_W1 <: ReproducingKernel_1

Define a type of reproducing kernel of Sobolev space $W^2_1 [0,1]$:

\[V(\eta, \xi) = \begin{cases} 1 + \eta \, , & 0 \le \eta \le \xi \le 1 \\ 1 + \xi \, , & 0 \le \xi \le \eta \le 1 \end{cases}\]
source
NormalSplines.RK_W2Type

struct RK_W2 <: ReproducingKernel_2

Define a type of reproducing kernel of Sobolev space $W^2_2 [0,1]$:

\[V(\eta, \xi) = \begin{cases} 1 + (\eta + \eta^2 / 2)\xi - \eta^3 / 6 \, , & 0 \le \eta \le \xi \le 1 \\ 1 + (\xi + \xi^2 / 2)\eta - \xi^3 / 6 \, , & 0 \le \xi \le \eta \le 1 \end{cases}\]
source
NormalSplines.RK_W3Type

struct RK_W3 <: ReproducingKernel_3

Define a type of reproducing kernel of Sobolev space $W^2_3 [0,1]$:

\[V(\eta, \xi) = \begin{cases} \sum_{i=0}^2 \frac{\xi^i}{!i} \left ( \frac{\eta^i}{!i} + (-1)^{i} \frac {\eta^{5 - i}}{(5 - i)!} \right ) \, , & 0 \le \eta \le \xi \le 1 \\ \sum_{i=0}^2 \frac{\eta^i}{!i} \left ( \frac{\xi^i}{!i} + (-1)^{i} \frac {\xi^{5 - i}}{(5 - i)!} \right ) \, , & 0 \le \xi \le \eta \le 1 \end{cases}\]
source

Bessel Potential space Reproducing Kernels

NormalSplines.RK_H1Type

struct RK_H1{T <: AbstractFloat} <: ReproducingKernel_1

Define a type of reproducing kernel of Bessel Potential space $H^1_ε (R)$:

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

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

struct RK_H2{T <: AbstractFloat} <: ReproducingKernel_2

Define a type of reproducing kernel of Bessel Potential space $H^2_ε (R)$:

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

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.
source
NormalSplines.RK_H3Type

struct RK_H3{T <: AbstractFloat} <: ReproducingKernel_3

Define a type of reproducing kernel of Bessel Potential space $H^3_ε (R)$:

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

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