Inclusive relation between B-spline spaces

Thm. Inclusive relation between B-spline spaces

For non-degenerate B-spline spaces, the following relationship holds.

\[\mathcal{P}[p,k] \subseteq \mathcal{P}[p',k'] \Leftrightarrow (m=p'-p \ge 0 \ \text{and} \ k+m\widehat{k}\subseteq k')\]

Base.issubsetMethod

Check inclusive relationship between B-spline spaces.

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

Examples

julia> P1 = BSplineSpace{1}(KnotVector([1,3,5,8]));


julia> P2 = BSplineSpace{1}(KnotVector([1,3,5,6,8,9]));


julia> P3 = BSplineSpace{2}(KnotVector([1,1,3,3,5,5,8,8]));


julia> P1 ⊆ P2
true

julia> P1 ⊆ P3
true

julia> P2 ⊆ P3
false

julia> P2 ⊈ P3
true
source

Here are plots of the B-spline basis functions of the spaces P1, P2, P3.

julia> P1 = BSplineSpace{1}(KnotVector([1,3,5,8]))BSplineSpace{1, Int64, KnotVector{Int64}}(KnotVector([1, 3, 5, 8]))
julia> P2 = BSplineSpace{1}(KnotVector([1,3,5,6,8,9]))BSplineSpace{1, Int64, KnotVector{Int64}}(KnotVector([1, 3, 5, 6, 8, 9]))
julia> P3 = BSplineSpace{2}(KnotVector([1,1,3,3,5,5,8,8]))BSplineSpace{2, Int64, KnotVector{Int64}}(KnotVector([1, 1, 3, 3, 5, 5, 8, 8]))
julia> plot( plot([t->bsplinebasis₊₀(P1,i,t) for i in 1:dim(P1)], 1, 9, ylims=(0,1), legend=false), plot([t->bsplinebasis₊₀(P2,i,t) for i in 1:dim(P2)], 1, 9, ylims=(0,1), legend=false), plot([t->bsplinebasis₊₀(P3,i,t) for i in 1:dim(P3)], 1, 9, ylims=(0,1), legend=false), layout=(3,1), link=:x )Plot{Plots.PlotlyBackend() n=11}

This means, there exists a $n \times n'$ matrix $A$ which holds:

\[\begin{aligned} B_{(i,p,k)} &=\sum_{j}A_{ij} B_{(j,p',k')} \\ n&=\dim(\mathcal{P}[p,k]) \\ n'&=\dim(\mathcal{P}[p',k']) \end{aligned}\]

You can calculate the change of basis matrix $A$ with changebasis.

julia> A12 = changebasis(P1,P2)2×4 SparseArrays.SparseMatrixCSC{Float64, Int32} with 3 stored entries:
 1.0   ⋅    ⋅         ⋅
  ⋅   1.0  0.666667   ⋅
julia> A13 = changebasis(P1,P3)2×5 SparseArrays.SparseMatrixCSC{Float64, Int32} with 6 stored entries: 0.5 1.0 0.5 ⋅ ⋅ ⋅ ⋅ 0.5 1.0 0.5
julia> plot(
           plot([t->bsplinebasis₊₀(P1,i,t) for i in 1:dim(P1)], 1, 9, ylims=(0,1), legend=false),
           plot([t->sum(A12[i,j]*bsplinebasis₊₀(P2,j,t) for j in 1:dim(P2)) for i in 1:dim(P1)], 1, 9, ylims=(0,1), legend=false),
           plot([t->sum(A13[i,j]*bsplinebasis₊₀(P3,j,t) for j in 1:dim(P3)) for i in 1:dim(P1)], 1, 9, ylims=(0,1), legend=false),
           layout=(3,1),
           link=:x
       )Plot{Plots.PlotlyBackend() n=6}
BasicBSpline.changebasis_RFunction

Return a coefficient matrix $A$ which satisfy

\[B_{(i,p,k)} = \sum_{j}A_{i,j}B_{(j,p',k')}\]

Assumption:

  • $P ⊆ P^{\prime}$
source
BasicBSpline.changebasis_IFunction

Return a coefficient matrix $A$ which satisfy

\[B_{(i,p,k)} = \sum_{j}A_{i,j}B_{(j,p',k')}\]

Assumption:

  • $P ⊑ P^{\prime}$
source
BasicBSpline.issqsubsetFunction

Check inclusive relationship between B-spline spaces.

\[\mathcal{P}[p,k] \sqsubseteq\mathcal{P}[p',k'] \Leftrightarrow \mathcal{P}[p,k]|_{[k_{p+1},k_{l-p}]} \subseteq\mathcal{P}[p',k']|_{[k'_{p'+1},k'_{l'-p'}]}\]

source
BasicBSpline.expandspaceFunction

Expand B-spline space with given additional degree and knotvector. The behavior of expandspace is same as expandspace_I.

source
BasicBSpline.expandspace_RFunction

Expand B-spline space with given additional degree and knotvector. This function is compatible with issubset ().

Examples

julia> k = KnotVector([0.0, 1.5, 2.5, 5.5, 8.0, 9.0, 9.5, 10.0]);


julia> P = BSplineSpace{2}(k);


julia> P′ = expandspace_R(P, Val(1), KnotVector([6.0]))
BSplineSpace{3, Float64, KnotVector{Float64}}(KnotVector([0.0, 0.0, 1.5, 1.5, 2.5, 2.5, 5.5, 5.5, 6.0, 8.0, 8.0, 9.0, 9.0, 9.5, 9.5, 10.0, 10.0]))

julia> P ⊆ P′
true

julia> P ⊑ P′
false

julia> domain(P)
2.5 .. 9.0

julia> domain(P′)
1.5 .. 9.5
source
BasicBSpline.expandspace_IFunction

Expand B-spline space with given additional degree and knotvector. This function is compatible with issqsubset ()

Examples

julia> k = KnotVector([0.0, 1.5, 2.5, 5.5, 8.0, 9.0, 9.5, 10.0]);


julia> P = BSplineSpace{2}(k);


julia> P′ = expandspace_I(P, Val(1), KnotVector([6.0]))
BSplineSpace{3, Float64, KnotVector{Float64}}(KnotVector([0.0, 1.5, 2.5, 2.5, 5.5, 5.5, 6.0, 8.0, 8.0, 9.0, 9.0, 9.5, 10.0]))

julia> P ⊆ P′
false

julia> P ⊑ P′
true

julia> domain(P)
2.5 .. 9.0

julia> domain(P′)
2.5 .. 9.0
source