B-spline space

Defnition

Before defining B-spline space, we'll define polynomial space with degree $p$.

Def. Polynomial space

Polynomial space with degree $p$.

\[\mathcal{P}[p] =\left\{f:\mathbb{R}\to\mathbb{R}\ ;\ t\mapsto a_0+a_1t^1+\cdots+a_pt^p \ \left| \ % a_i\in \mathbb{R} \right. \right\}\]

This space $\mathcal{P}[p]$ is a $(p+1)$-dimensional linear space.

Note that $\{t\mapsto t^i\}_{0 \le i \le p}$ is a basis of $\mathcal{P}[p]$, and also the set of Bernstein polynomial $\{B_{(i,p)}\}_i$ is a basis of $\mathcal{P}[p]$.

\[\begin{aligned} B_{(i,p)}(t) &=\binom{p}{i-1}t^{i-1}(1-t)^{p-i+1} &(i=1, \dots, p+1) \end{aligned}\]

Where $\binom{p}{i-1}$ is a binomial coefficient. You can try Bernstein polynomial on desmos graphing calculator!

Def. B-spline space

For given polynomial degree $p\ge 0$ and knot vector $k=(k_1,\dots,k_l)$, B-spline space $\mathcal{P}[p,k]$ is defined as follows:

\[\mathcal{P}[p,k] =\left\{f:\mathbb{R}\to\mathbb{R} \ \left| \ % \begin{gathered} \operatorname{supp}(f)\subseteq [k_1, k_l] \\ \exists \tilde{f}\in\mathcal{P}[p], f|_{[k_{i}, k_{i+1})} = \tilde{f}|_{[k_{i}, k_{i+1})} \\ \forall t \in \mathbb{R}, \exists \delta > 0, f|_{(t-\delta,t+\delta)}\in C^{p-\mathfrak{n}_k(t)} \end{gathered} \right. \right\}\]

Note that each element of the space $\mathcal{P}[p,k]$ is a piecewise polynomial.

[TODO: fig]

BasicBSpline.BSplineSpaceType

Construct B-spline space from given polynominal degree and knot vector.

\[\mathcal{P}[p,k]\]

Examples

julia> p = 2
2

julia> k = KnotVector([1,3,5,6,8,9])
KnotVector([1, 3, 5, 6, 8, 9])

julia> BSplineSpace{p}(k)
BSplineSpace{2, Int64, KnotVector{Int64}}(KnotVector([1, 3, 5, 6, 8, 9]))
source

Degeneration

Def. Degeneration

A B-spline space is said to be non-degenerate if its degree and knot vector satisfies following property:

\[\begin{aligned} k_{i}&<k_{i+p+1} & (1 \le i \le l-p-1) \end{aligned}\]

BasicBSpline.isnondegenerateFunction

Check if given B-spline space is non-degenerate.

Examples

julia> isnondegenerate(BSplineSpace{2}(KnotVector([1,3,5,6,8,9])))
true

julia> isnondegenerate(BSplineSpace{1}(KnotVector([1,3,3,3,8,9])))
false
source
BasicBSpline.isdegenerateMethod

Check if given B-spline space is degenerate.

Examples

julia> isdegenerate(BSplineSpace{2}(KnotVector([1,3,5,6,8,9])))
false

julia> isdegenerate(BSplineSpace{1}(KnotVector([1,3,3,3,8,9])))
true
source

Dimensions

Thm. Dimension of B-spline space

The B-spline space is a linear space, and if a B-spline space is non-degenerate, its dimension is calculated by:

\[\dim(\mathcal{P}[p,k])=\# k - p -1\]

BasicBSpline.dimFunction

Return dimention of a B-spline space.

\[\dim(\mathcal{P}[p,k]) =\# k - p -1\]

Examples

julia> dim(BSplineSpace{1}(KnotVector([1,2,3,4,5,6,7])))
5

julia> dim(BSplineSpace{1}(KnotVector([1,2,4,4,4,6,7])))
5

julia> dim(BSplineSpace{1}(KnotVector([1,2,3,5,5,5,7])))
5
source
BasicBSpline.exactdim_RMethod

Exact dimension of a B-spline space.

Examples

julia> exactdim_R(BSplineSpace{1}(KnotVector([1,2,3,4,5,6,7])))
5

julia> exactdim_R(BSplineSpace{1}(KnotVector([1,2,4,4,4,6,7])))
4

julia> exactdim_R(BSplineSpace{1}(KnotVector([1,2,3,5,5,5,7])))
4
source
BasicBSpline.exactdim_IMethod

Exact dimension of a B-spline space.

Examples

julia> exactdim_I(BSplineSpace{1}(KnotVector([1,2,3,4,5,6,7])))
5

julia> exactdim_I(BSplineSpace{1}(KnotVector([1,2,4,4,4,6,7])))
4

julia> exactdim_I(BSplineSpace{1}(KnotVector([1,2,3,5,5,5,7])))
3
source