API

BasicBSpline.jl

BasicBSpline.BSplineDerivativeSpaceType
BSplineDerivativeSpace{r}(P::BSplineSpace)

Construct derivative of B-spline space from given differential order and B-spline space.

\[D^{r}(\mathcal{P}[p,k]) =\left\{t \mapsto \left. \frac{d^r f}{dt^r}(t) \ \right| \ f \in \mathcal{P}[p,k] \right\}\]

Examples

julia> P = BSplineSpace{2}(KnotVector([1,2,3,4,5,6]))
BSplineSpace{2, Int64, KnotVector{Int64}}(KnotVector([1, 2, 3, 4, 5, 6]))

julia> dP = BSplineDerivativeSpace{1}(P)
BSplineDerivativeSpace{1, BSplineSpace{2, Int64, KnotVector{Int64}}, Int64}(BSplineSpace{2, Int64, KnotVector{Int64}}(KnotVector([1, 2, 3, 4, 5, 6])))

julia> degree(P), degree(dP)
(2, 1)
source
BasicBSpline.BSplineManifoldType

Construct B-spline manifold from given control points and B-spline spaces.

Examples

julia> using StaticArrays

julia> P = BSplineSpace{2}(KnotVector([0,0,0,1,1,1]))
BSplineSpace{2, Int64, KnotVector{Int64}}(KnotVector([0, 0, 0, 1, 1, 1]))

julia> a = [SVector(1,0), SVector(1,1), SVector(0,1)]
3-element Vector{SVector{2, Int64}}:
 [1, 0]
 [1, 1]
 [0, 1]

julia> M = BSplineManifold(a, P);


julia> M(0.4)
2-element SVector{2, Float64} with indices SOneTo(2):
 0.84
 0.64

julia> M(1.2)
ERROR: DomainError with 1.2:
The input 1.2 is out of domain 0 .. 1.
[...]
source
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
BasicBSpline.EmptyKnotVectorType

Knot vector with zero-element.

\[k=()\]

This struct is intended for internal use.

Examples

julia> EmptyKnotVector()
EmptyKnotVector{Bool}()

julia> EmptyKnotVector{Float64}()
EmptyKnotVector{Float64}()
source
BasicBSpline.KnotVectorType

Construct knot vector from given array.

\[k=(k_1,\dots,k_l)\]

Examples

julia> k = KnotVector([1,2,3])
KnotVector([1, 2, 3])

julia> k = KnotVector(1:3)
KnotVector([1, 2, 3])
source
BasicBSpline.RationalBSplineManifoldType

Construct Rational B-spline manifold from given control points, weights and B-spline spaces.

Examples

julia> using StaticArrays, LinearAlgebra

julia> P = BSplineSpace{2}(KnotVector([0,0,0,1,1,1]))
BSplineSpace{2, Int64, KnotVector{Int64}}(KnotVector([0, 0, 0, 1, 1, 1]))

julia> w = [1, 1/√2, 1]
3-element Vector{Float64}:
 1.0
 0.7071067811865475
 1.0

julia> a = [SVector(1,0), SVector(1,1), SVector(0,1)]
3-element Vector{SVector{2, Int64}}:
 [1, 0]
 [1, 1]
 [0, 1]

julia> M = RationalBSplineManifold(a,w,P);  # 1/4 arc


julia> M(0.3)
2-element SVector{2, Float64} with indices SOneTo(2):
 0.8973756499953727
 0.4412674277525845

julia> norm(M(0.3))
1.0
source
BasicBSpline.SubKnotVectorType

A type to represetnt sub knot vector.

\[k=(k_1,\dots,k_l)\]

Examples

julia> k = knotvector"1 11 211"
KnotVector([1, 3, 4, 6, 6, 7, 8])

julia> view(k, 2:5)
SubKnotVector([3, 4, 6, 6])
source
BasicBSpline.UniformKnotVectorType

Construct uniform knot vector from given range.

\[k=(k_1,\dots,k_l)\]

Examples

julia> k = UniformKnotVector(1:8)
UniformKnotVector(1:8)

julia> UniformKnotVector(8:-1:3)
UniformKnotVector(3:1:8)
source
Base.:*Method

Product of integer and knot vector

\[\begin{aligned} m\cdot k&=\underbrace{k+\cdots+k}_{m} \end{aligned}\]

For example, $2\cdot (1,2,2,5)=(1,1,2,2,2,2,5,5)$.

Examples

julia> k = KnotVector([1,2,2,5]);


julia> 2 * k
KnotVector([1, 1, 2, 2, 2, 2, 5, 5])
source
Base.:+Method

Sum of knot vectors

\[\begin{aligned} k^{(1)}+k^{(2)} &=(k^{(1)}_1, \dots, k^{(1)}_{l^{(1)}}) + (k^{(2)}_1, \dots, k^{(2)}_{l^{(2)}}) \\ &=(\text{sort of union of} \ k^{(1)} \ \text{and} \ k^{(2)} \text{)} \end{aligned}\]

For example, $(1,2,3,5)+(4,5,8)=(1,2,3,4,5,5,8)$.

Examples

julia> k1 = KnotVector([1,2,3,5]);


julia> k2 = KnotVector([4,5,8]);


julia> k1 + k2
KnotVector([1, 2, 3, 4, 5, 5, 8])
source
Base.issubsetMethod

Check a inclusive relationship $k\subseteq k'$, for example:

\[\begin{aligned} (1,2) &\subseteq (1,2,3) \\ (1,2,2) &\not\subseteq (1,2,3) \\ (1,2,3) &\subseteq (1,2,3) \\ \end{aligned}\]

Examples

julia> KnotVector([1,2]) ⊆ KnotVector([1,2,3])
true

julia> KnotVector([1,2,2]) ⊆ KnotVector([1,2,3])
false

julia> KnotVector([1,2,3]) ⊆ KnotVector([1,2,3])
true
source
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
Base.lengthMethod

Length of knot vector

\[\begin{aligned} \#{k} &=(\text{number of knot elements of} \ k) \\ \end{aligned}\]

For example, $\#{(1,2,2,3)}=4$.

Examples

julia> k = KnotVector([1,2,2,3]);


julia> length(k)
4
source
Base.uniqueMethod

Unique elements of knot vector.

\[\begin{aligned} \widehat{k} &=(\text{unique knot elements of} \ k) \\ \end{aligned}\]

For example, $\widehat{(1,2,2,3)}=(1,2,3)$.

Examples

julia> k = KnotVector([1,2,2,3]);


julia> unique(k)
KnotVector([1, 2, 3])
source
BasicBSpline.__changebasis_IFunction
__changebasis_I(P::AbstractFunctionSpace, P′::AbstractFunctionSpace)

Internal function for changebasis_I.

Implicit assumption:

  • P ⊑ P′
  • isnondegenerate_I(P′, 1)
  • isnondegenerate_I(P′, dim(P′))
source
BasicBSpline._lower_RFunction

Internal methods for obtaining a B-spline space with one degree lower.

\[\begin{aligned} \mathcal{P}[p,k] &\mapsto \mathcal{P}[p-1,k] \\ D^r\mathcal{P}[p,k] &\mapsto D^{r-1}\mathcal{P}[p-1,k] \end{aligned}\]

source
BasicBSpline.bsplinebasisMethod

$i$-th B-spline basis function. Modified version.

\[\begin{aligned} {B}_{(i,p,k)}(t) &= \frac{t-k_{i}}{k_{i+p}-k_{i}}{B}_{(i,p-1,k)}(t) +\frac{k_{i+p+1}-t}{k_{i+p+1}-k_{i+1}}{B}_{(i+1,p-1,k)}(t) \\ {B}_{(i,0,k)}(t) &= \begin{cases} &1\quad (k_{i} \le t < k_{i+1}) \\ &1\quad (k_{i} < t = k_{i+1}=k_{l}) \\ &0\quad (\text{otherwise}) \end{cases} \end{aligned}\]

Examples

julia> P = BSplineSpace{0}(KnotVector(1:6))
BSplineSpace{0, Int64, KnotVector{Int64}}(KnotVector([1, 2, 3, 4, 5, 6]))

julia> bsplinebasis.(P,1:5,(1:6)')
5×6 Matrix{Float64}:
 1.0  0.0  0.0  0.0  0.0  0.0
 0.0  1.0  0.0  0.0  0.0  0.0
 0.0  0.0  1.0  0.0  0.0  0.0
 0.0  0.0  0.0  1.0  0.0  0.0
 0.0  0.0  0.0  0.0  1.0  1.0
source
BasicBSpline.bsplinebasisallFunction

B-spline basis functions at point t on i-th interval.

Examples

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

julia> p = 2
2

julia> P = BSplineSpace{p}(k)
BSplineSpace{2, Float64, KnotVector{Float64}}(KnotVector([0.0, 1.5, 2.5, 5.5, 8.0, 9.0, 9.5, 10.0]))

julia> t = 5.7
5.7

julia> i = intervalindex(P,t)
2

julia> bsplinebasisall(P,i,t)
3-element SVector{3, Float64} with indices SOneTo(3):
 0.3847272727272727
 0.6107012987012989
 0.00457142857142858

julia> bsplinebasis.(P,i:i+p,t)
3-element Vector{Float64}:
 0.38472727272727264
 0.6107012987012989
 0.00457142857142858
source
BasicBSpline.bsplinebasis′Function
bsplinebasis′(::AbstractFunctionSpace, ::Integer, ::Real) -> Real

1st derivative of B-spline basis function. Modified version.

\[\dot{B}_{(i,p,k)}(t) =p\left(\frac{1}{k_{i+p}-k_{i}}B_{(i,p-1,k)}(t)-\frac{1}{k_{i+p+1}-k_{i+1}}B_{(i+1,p-1,k)}(t)\right)\]

bsplinebasis′(P, i, t) is equivalent to bsplinebasis(derivative(P), i, t).

source
BasicBSpline.bsplinebasis′₊₀Function
bsplinebasis′₊₀(::AbstractFunctionSpace, ::Integer, ::Real) -> Real

1st derivative of B-spline basis function. Right-sided limit version.

\[\dot{B}_{(i,p,k)}(t) =p\left(\frac{1}{k_{i+p}-k_{i}}B_{(i,p-1,k)}(t)-\frac{1}{k_{i+p+1}-k_{i+1}}B_{(i+1,p-1,k)}(t)\right)\]

bsplinebasis′₊₀(P, i, t) is equivalent to bsplinebasis₊₀(derivative(P), i, t).

source
BasicBSpline.bsplinebasis′₋₀Function
bsplinebasis′₋₀(::AbstractFunctionSpace, ::Integer, ::Real) -> Real

1st derivative of B-spline basis function. Left-sided limit version.

\[\dot{B}_{(i,p,k)}(t) =p\left(\frac{1}{k_{i+p}-k_{i}}B_{(i,p-1,k)}(t)-\frac{1}{k_{i+p+1}-k_{i+1}}B_{(i+1,p-1,k)}(t)\right)\]

bsplinebasis′₋₀(P, i, t) is equivalent to bsplinebasis₋₀(derivative(P), i, t).

source
BasicBSpline.bsplinebasis₊₀Method

$i$-th B-spline basis function. Right-sided limit version.

\[\begin{aligned} {B}_{(i,p,k)}(t) &= \frac{t-k_{i}}{k_{i+p}-k_{i}}{B}_{(i,p-1,k)}(t) +\frac{k_{i+p+1}-t}{k_{i+p+1}-k_{i+1}}{B}_{(i+1,p-1,k)}(t) \\ {B}_{(i,0,k)}(t) &= \begin{cases} &1\quad (k_{i}\le t< k_{i+1})\\ &0\quad (\text{otherwise}) \end{cases} \end{aligned}\]

Examples

julia> P = BSplineSpace{0}(KnotVector(1:6))
BSplineSpace{0, Int64, KnotVector{Int64}}(KnotVector([1, 2, 3, 4, 5, 6]))

julia> bsplinebasis₊₀.(P,1:5,(1:6)')
5×6 Matrix{Float64}:
 1.0  0.0  0.0  0.0  0.0  0.0
 0.0  1.0  0.0  0.0  0.0  0.0
 0.0  0.0  1.0  0.0  0.0  0.0
 0.0  0.0  0.0  1.0  0.0  0.0
 0.0  0.0  0.0  0.0  1.0  0.0
source
BasicBSpline.bsplinebasis₋₀Method

$i$-th B-spline basis function. Left-sided limit version.

\[\begin{aligned} {B}_{(i,p,k)}(t) &= \frac{t-k_{i}}{k_{i+p}-k_{i}}{B}_{(i,p-1,k)}(t) +\frac{k_{i+p+1}-t}{k_{i+p+1}-k_{i+1}}{B}_{(i+1,p-1,k)}(t) \\ {B}_{(i,0,k)}(t) &= \begin{cases} &1\quad (k_{i}< t\le k_{i+1})\\ &0\quad (\text{otherwise}) \end{cases} \end{aligned}\]

Examples

julia> P = BSplineSpace{0}(KnotVector(1:6))
BSplineSpace{0, Int64, KnotVector{Int64}}(KnotVector([1, 2, 3, 4, 5, 6]))

julia> bsplinebasis₋₀.(P,1:5,(1:6)')
5×6 Matrix{Float64}:
 0.0  1.0  0.0  0.0  0.0  0.0
 0.0  0.0  1.0  0.0  0.0  0.0
 0.0  0.0  0.0  1.0  0.0  0.0
 0.0  0.0  0.0  0.0  1.0  0.0
 0.0  0.0  0.0  0.0  0.0  1.0
source
BasicBSpline.bsplinebasis₋₀IMethod

$i$-th B-spline basis function. Modified version (2).

\[\begin{aligned} {B}_{(i,p,k)}(t) &= \frac{t-k_{i}}{k_{i+p}-k_{i}}{B}_{(i,p-1,k)}(t) +\frac{k_{i+p+1}-t}{k_{i+p+1}-k_{i+1}}{B}_{(i+1,p-1,k)}(t) \\ {B}_{(i,0,k)}(t) &= \begin{cases} &1\quad (k_{1+p} \le k_{i} < t \le k_{i+1}) \\ &1\quad (t = k_{1+p} = k_{i} < k_{i+1}) \\ &1\quad (k_{i} \le t < k_{i+1} \le k_{1+p}) \\ &0\quad (\text{otherwise}) \end{cases} \end{aligned}\]

Examples

julia> P = BSplineSpace{0}(KnotVector(1:6))
BSplineSpace{0, Int64, KnotVector{Int64}}(KnotVector([1, 2, 3, 4, 5, 6]))

julia> BasicBSpline.bsplinebasis₋₀I.(P,1:5,(1:6)')
5×6 Matrix{Float64}:
 1.0  1.0  0.0  0.0  0.0  0.0
 0.0  0.0  1.0  0.0  0.0  0.0
 0.0  0.0  0.0  1.0  0.0  0.0
 0.0  0.0  0.0  0.0  1.0  0.0
 0.0  0.0  0.0  0.0  0.0  1.0
source
BasicBSpline.bsplinesupportFunction

Return the support of $i$-th B-spline basis function.

\[\operatorname{supp}(B_{(i,p,k)})=[k_{i},k_{i+p+1}]\]

Examples

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

julia> P = BSplineSpace{2}(k)
BSplineSpace{2, Float64, KnotVector{Float64}}(KnotVector([0.0, 1.5, 2.5, 5.5, 8.0, 9.0, 9.5, 10.0]))

julia> bsplinesupport(P,1)
0.0 .. 5.5

julia> bsplinesupport(P,2)
1.5 .. 8.0
source
BasicBSpline.bsplinesupport_IMethod

Return the support of $i$-th B-spline basis function.

\[\operatorname{supp}(B_{(i,p,k)}|_I)=[k_{i},k_{i+p+1}] \cap I \qquad (I = [k_{1+p}, k_{l-p}])\]

Examples

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

julia> P = BSplineSpace{2}(k)
BSplineSpace{2, Float64, KnotVector{Float64}}(KnotVector([0.0, 1.5, 2.5, 5.5, 8.0, 9.0, 9.5, 10.0]))

julia> bsplinesupport_I(P,1)
2.5 .. 5.5

julia> bsplinesupport_I(P,2)
2.5 .. 8.0
source
BasicBSpline.bsplinesupport_RMethod

Return the support of $i$-th B-spline basis function.

\[\operatorname{supp}(B_{(i,p,k)})=[k_{i},k_{i+p+1}]\]

Examples

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

julia> P = BSplineSpace{2}(k)
BSplineSpace{2, Float64, KnotVector{Float64}}(KnotVector([0.0, 1.5, 2.5, 5.5, 8.0, 9.0, 9.5, 10.0]))

julia> bsplinesupport_R(P,1)
0.0 .. 5.5

julia> bsplinesupport_R(P,2)
1.5 .. 8.0
source
BasicBSpline.changebasisMethod
changebasis(P::AbstractFunctionSpace, P′::AbstractFunctionSpace)

Return changebasis_R(P, P′) if $P ⊆ P′$, otherwise changebasis_R(P, P′) if $P ⊑ P′$. Throw an error if $P ⊈ P′$ and $P ⋢ P′$.

source
BasicBSpline.changebasis_IMethod
changebasis_I(P::AbstractFunctionSpace, P′::AbstractFunctionSpace)

Return a coefficient matrix $A$ which satisfy

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

Examples

julia> P = BSplineSpace{2}(knotvector"3 13")
BSplineSpace{2, Int64, KnotVector{Int64}}(KnotVector([1, 1, 1, 3, 4, 4, 4]))

julia> P′ = BSplineSpace{3}(knotvector"4124")
BSplineSpace{3, Int64, KnotVector{Int64}}(KnotVector([1, 1, 1, 1, 2, 3, 3, 4, 4, 4, 4]))

julia> P ⊑ P′
true

julia> changebasis_I(P, P′)
4×7 SparseArrays.SparseMatrixCSC{Float64, Int32} with 13 stored entries:
 1.0  0.666667  0.166667   ⋅         ⋅         ⋅         ⋅
  ⋅   0.333333  0.722222  0.555556  0.111111   ⋅         ⋅
  ⋅    ⋅        0.111111  0.444444  0.888889  0.666667   ⋅
  ⋅    ⋅         ⋅         ⋅         ⋅        0.333333  1.0
source
BasicBSpline.changebasis_RMethod
changebasis_R(P::AbstractFunctionSpace, P′::AbstractFunctionSpace)

Return a coefficient matrix $A$ which satisfy

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

Examples

julia> P = BSplineSpace{2}(knotvector"3 13")
BSplineSpace{2, Int64, KnotVector{Int64}}(KnotVector([1, 1, 1, 3, 4, 4, 4]))

julia> P′ = BSplineSpace{3}(knotvector"4124")
BSplineSpace{3, Int64, KnotVector{Int64}}(KnotVector([1, 1, 1, 1, 2, 3, 3, 4, 4, 4, 4]))

julia> P ⊆ P′
true

julia> changebasis_R(P, P′)
4×7 SparseArrays.SparseMatrixCSC{Float64, Int32} with 13 stored entries:
 1.0  0.666667  0.166667   ⋅         ⋅         ⋅         ⋅
  ⋅   0.333333  0.722222  0.555556  0.111111   ⋅         ⋅
  ⋅    ⋅        0.111111  0.444444  0.888889  0.666667   ⋅
  ⋅    ⋅         ⋅         ⋅         ⋅        0.333333  1.0
source
BasicBSpline.countknotsMethod

For given knot vector $k$, the following function $\mathfrak{n}_k:\mathbb{R}\to\mathbb{Z}$ represents the number of knots that duplicate the knot vector $k$.

\[\mathfrak{n}_k(t) = \#\{i \mid k_i=t \}\]

For example, if $k=(1,2,2,3)$, then $\mathfrak{n}_k(0.3)=0$, $\mathfrak{n}_k(1)=1$, $\mathfrak{n}_k(2)=2$.

julia> k = KnotVector([1,2,2,3]);


julia> countknots(k,0.3)
0

julia> countknots(k,1.0)
1

julia> countknots(k,2.0)
2
source
BasicBSpline.derivativeMethod
derivative(::BSplineDerivativeSpace{r}) -> BSplineDerivativeSpace{r+1}
derivative(::BSplineSpace) -> BSplineDerivativeSpace{1}

Derivative of B-spline related space.

Examples

julia> BSplineSpace{2}(KnotVector(0:5))
BSplineSpace{2, Int64, KnotVector{Int64}}(KnotVector([0, 1, 2, 3, 4, 5]))

julia> BasicBSpline.derivative(ans)
BSplineDerivativeSpace{1, BSplineSpace{2, Int64, KnotVector{Int64}}, Int64}(BSplineSpace{2, Int64, KnotVector{Int64}}(KnotVector([0, 1, 2, 3, 4, 5])))

julia> BasicBSpline.derivative(ans)
BSplineDerivativeSpace{2, BSplineSpace{2, Int64, KnotVector{Int64}}, Int64}(BSplineSpace{2, Int64, KnotVector{Int64}}(KnotVector([0, 1, 2, 3, 4, 5])))
source
BasicBSpline.dimMethod

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_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
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.expandspaceMethod

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

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

Return an index of a interval in the domain of B-spline space

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> domain(P)
2.5 .. 9.0

julia> intervalindex(P,2.6)
1

julia> intervalindex(P,5.6)
2

julia> intervalindex(P,8.5)
3

julia> intervalindex(P,9.5)
3
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
BasicBSpline.isnondegenerateMethod

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

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

Calculate $r$-nomial coefficient

r_nomial(n, k, r)

This function is considered as internal.

\[(1+x+\cdots+x^r)^n = \sum_{k} a_{n,k,r} x^k\]

source
BasicBSpline.unbounded_mappingFunction
unbounded_mapping(M::BSplineManifold{Dim}, t::Vararg{Real,Dim})

Examples

julia> P = BSplineSpace{1}(KnotVector([0,0,1,1]))
BSplineSpace{1, Int64, KnotVector{Int64}}(KnotVector([0, 0, 1, 1]))

julia> domain(P)
0 .. 1

julia> M = BSplineManifold([0,1], P);


julia> unbounded_mapping(M, 0.1)
0.1

julia> M(0.1)
0.1

julia> unbounded_mapping(M, 1.2)
1.2

julia> M(1.2)
ERROR: DomainError with 1.2:
The input 1.2 is out of domain 0 .. 1.
[...]
source
BasicBSpline.@knotvector_strMacro
@knotvector_str -> KnotVector

Construct a knotvector by specifying the numbers of duplicates of knots.

Examples

julia> knotvector"11111"
KnotVector([1, 2, 3, 4, 5])

julia> knotvector"123"
KnotVector([1, 2, 2, 3, 3, 3])

julia> knotvector" 2 2 2"
KnotVector([2, 2, 4, 4, 6, 6])

julia> knotvector"     1"
KnotVector([6])
source

BasicBSplineFitting.jl

BasicBSplineFitting.fittingcontrolpointsMethod

Fitting controlpoints with least squares method.

fittingcontrolpoints(func, Ps::Tuple)

This function will calculate $\bm{a}_i$ to minimize the following integral.

\[\int_I \left\|f(t)-\sum_i B_{(i,p,k)}(t) \bm{a}_i\right\|^2 dt\]

Similarly, for the two-dimensional case, minimize the following integral.

\[\int_{I^1 \times I^2} \left\|f(t^1, t^2)-\sum_{i,j} B_{(i,p^1,k^1)}(t^1)B_{(j,p^2,k^2)}(t^2) \bm{a}_{ij}\right\|^2 dt^1dt^2\]

Currently, this function supports up to three dimensions.

Examples

julia> f(t) = SVector(cos(t),sin(t),t);

julia> P = BSplineSpace{3}(KnotVector(range(0,2π,30)) + 3*KnotVector([0,2π]));

julia> a = fittingcontrolpoints(f, P);

julia> M = BSplineManifold(a, P);

julia> norm(M(1) - f(1)) < 1e-5
true
source