Knot sets
Three built-in knot sets are provided out-of-the-box:
Arbitrary knot set
Arbitrary knot sets allow custom placement of the individual knots, as well as custom multiplicity of interior knots. This can be used for problems where discontinuity at a certain location is desired.
CompactBases.ArbitraryKnotSet
— TypeArbitraryKnotSet{k,ml,mr}(t)
An arbitrary knot set of order k
and left and right multiplicities of ml
and mr
, respectively. The knot set is specified by the non-decreasing vector t
; the same knot may appear multiple times, which influences the continuity of the B-splines at that location.
CompactBases.ArbitraryKnotSet
— TypeArbitraryKnotset(k, t[, ml=k, mr=k])
Construct an order-k
arbitrary knot set, with the locations of the knots given by non-decreasing vector t
.
Linear knot set
CompactBases.LinearKnotSet
— TypeLinearKnotSet{k,ml,mr}(t)
A knot set of order k
and left and right multiplicities of ml
and mr
, respectively, whose knots are uniformly distributed according to the range t
. The interior basis functions are thus Cᵏ⁻²-continuous.
CompactBases.LinearKnotSet
— TypeLinearKnotSet(k, a, b, N[, ml=k, mr=k])
Construct an order-k
linear knot set spanning from a
to b
, with N
intervals.
Examples
julia> LinearKnotSet(2, 0, 1, 3)
6-element LinearKnotSet{2,2,2,Float64,StepRangeLen{Float64,Base.TwicePrecision{Float64},Base.TwicePrecision{Float64}}}:
0.0
0.0
0.3333333333333333
0.6666666666666666
1.0
1.0
Exponential knot set
Exponential knot sets are useful for approximating exponentially varying functions (e.g. bound states of atoms). The quadrature points of each interval are distributed as were the knot set piecewise linear.
CompactBases.ExpKnotSet
— TypeExpKnotSet{k,ml,mr}(exponents, base, t, include0)
A knot set of order k
and left and right multiplicities of ml
and mr
, respectively, whose knots are exponentially distributed according to t = base .^ exponents
, optionally including 0 as the left endpoint.
CompactBases.ExpKnotSet
— MethodExpKnotSet(k, a, b, N[, ml=k, mr=k, base=10, include0=true])
Construct an order-k
knot spanning from base^a
to base^b
in N
intervals, optionally including 0 as the left endpoint.
Examples
julia> ExpKnotSet(2, -4, 2, 7)
10-element ExpKnotSet{2,2,2,Float64,StepRangeLen{Float64,Base.TwicePrecision{Float64},Base.TwicePrecision{Float64}},Array{Float64,1}}:
0.0
0.0
0.0001
0.001
0.01
0.1
1.0
10.0
100.0
100.0
Reference
CompactBases.order
— Functionorder(t)
Returns the order k
of the knot set t
.
CompactBases.numintervals
— Functionnuminterval(t)
Returns the number of intervals generated by the knot set t
.
Examples
julia> numintervals(LinearKnotSet(3, 0, 1, 2))
2
CompactBases.numfunctions
— Functionnumfunctions(t)
Returns the number of basis functions generated by knot set t
.
Examples
julia> numfunctions(LinearKnotSet(3, 0, 1, 2))
4
julia> numfunctions(LinearKnotSet(5, 0, 1, 2))
6
Base.first
— Functionfirst(t)
Return the first knot of t
.
Base.last
— Functionlast(t)
Return the last knot of t
.
Base.length
— Functionlength(t)
Return the number of knots of t
.
Examples
julia> length(LinearKnotSet(3, 0, 1, 3))
8
julia> length(LinearKnotSet(3, 0, 1, 3, 1, 1))
4
Base.getindex
— Functiongetindex(t, i)
Return the i
th knot of the knot set t
, accounting for the multiplicities of the endpoints.
Examples
julia> LinearKnotSet(3, 0, 1, 3)[2]
0.0
julia> LinearKnotSet(3, 0, 1, 3, 1, 1)[2]
0.3333333333333333
CompactBases.nonempty_intervals
— Functionnonempty_intervals(t)
Return the indices of all intervals of the knot set t
that are non-empty.
Examples
julia> nonempty_intervals(ArbitraryKnotSet(3, [0.0, 1, 1, 3, 4, 6], 1, 3))
4-element Array{Int64,1}:
1
3
4
5
Internals
CompactBases.AbstractKnotSet
— TypeAbstractKnotSet{k,ml,mr,T}
Abstract base for B-spline knot sets. T
is the eltype
of the knot set, k
is the order of the piecewise polynomials (order = degree + 1
) and ml
and mr
are the knot multiplicities of the left and right endpoints, respectively.
CompactBases.assert_multiplicities
— Functionassert_multiplicities(k,ml,mr,t)
Assert that the multiplicities at the endpoints, ml
and mr
, respectively, are consistent with the order k
. Also check that the amount of knots in the knot set t
are enough to support the requested order k
.
CompactBases.leftmultiplicity
— Functionleftmultiplicity(t)
Return the multiplicity of the left endpoint.
CompactBases.rightmultiplicity
— Functionrightmultiplicity(t)
Return the multiplicity of the right endpoint.
CompactBases.find_interval
— Functionfind_interval(t, x[, i=ml])
Find the interval in the knot set t
that includes x
, starting from interval i
(which by default is the first non-zero interval of the knot set). The search complexity is linear, but by storing the result and using it as starting point for the next call to find_interval
, the knot set need only be traversed once.
CompactBases.within_interval
— Functionwithin_interval(x, interval)
Return the indices of the elements of x
that lie within the given interval
.
CompactBases.within_support
— Functionwithin_support(x, t, j)
Return the indices of the elements of x
that lie withing the compact support of the j
th basis function (enumerated 1..n
), given the knot set t
. For each index of x
that is covered, the index k
of the interval within which x[i]
falls is also returned.
Quadrature functions
Gauß–Legendre quadrature can be used to exactly calculate integrals of polynomials, and very accurately the integrals of smoothly varying functions. The approximation is given by
\[\begin{equation} \int\limits_a^b \diff{x}f(x)\approx \frac{b-a}{2}\sum_{i=1}^n w_i f\left(\frac{b-a}{2}x_i+\frac{a+b}{2}\right), \end{equation}\]
where $x_i$ are the roots of the quadrature and $w_i$ the corresponding weights, given on the elementary interval $[-1,1]$.
CompactBases.num_quadrature_points
— Functionnum_quadrature_points(k, k′)
The number of quadrature points needed to exactly compute the matrix elements of an operator of polynomial order k′
with respect to a basis of order k
.
Missing docstring for CompactBases.lgwt!
. Check Documenter's build log for details.
CompactBases.lgwt
— Functionlgwt(t, N) -> (x,w)
Generate the N
Gauß–Legendre quadrature roots x
and associated weights w
, with respect to the B-spline basis generated by the knot set t
.
Examples
julia> CompactBases.lgwt(LinearKnotSet(2, 0, 1, 3), 2)
([0.0704416, 0.262892, 0.403775, 0.596225, 0.737108, 0.929558], [0.166667, 0.166667, 0.166667, 0.166667, 0.166667, 0.166667])
julia> CompactBases.lgwt(ExpKnotSet(2, -4, 2, 7), 2)
([2.11325e-5, 7.88675e-5, 0.000290192, 0.000809808, 0.00290192, 0.00809808, 0.0290192, 0.0809808, 0.290192, 0.809808, 2.90192, 8.09808, 29.0192, 80.9808], [5.0e-5, 5.0e-5, 0.00045, 0.00045, 0.0045, 0.0045, 0.045, 0.045, 0.45, 0.45, 4.5, 4.5, 45.0, 45.0])