B-spline manifold
Multi-dimensional B-spline
The tensor product of B-spline spaces $\mathcal{P}[p^1,k^1]\otimes\mathcal{P}[p^2,k^2]$ is a linear space with the following basis.
\[\mathcal{P}[p^1,k^1]\otimes\mathcal{P}[p^2,k^2] = \operatorname*{span}_{i,j} (B_{(i,p^1,k^1)} \otimes B_{(j,p^2,k^2)})\]
where the basis are defined as
\[(B_{(i,p^1,k^1)} \otimes B_{(j,p^2,k^2)})(t^1, t^2) = B_{(i,p^1,k^1)}(t^1) \cdot B_{(j,p^2,k^2)}(t^2)\]
Higher dimensional tensor products $\mathcal{P}[p^1,k^1]\otimes\cdots\otimes\mathcal{P}[p^d,k^d]$ are defined similarly.
B-spline manifold
B-spline manifold is a parametric representation of a shape.
For given $d$-dimensional B-spline basis functions $B_{(i^1,p^1,k^1)} \otimes \cdots \otimes B_{(i^d,p^d,k^d)}$ and given points $\bm{a}_{i^1 \dots i^d} \in V$, B-spline manifold is defined by the following equality:
\[\bm{p}(t^1,\dots,t^d;\bm{a}_{i^1 \dots i^d}) =\sum_{i^1,\dots,i^d}(B_{(i^1,p^1,k^1)} \otimes \cdots \otimes B_{(i^d,p^d,k^d)})(t^1,\dots,t^d) \bm{a}_{i^1 \dots i^d}\]
Where $\bm{a}_{i^1 \dots i^d}$ are called control points.
We will also write $\bm{p}(t^1,\dots,t^d; \bm{a})$, $\bm{p}(t^1,\dots,t^d)$, $\bm{p}(t; \bm{a})$ or $\bm{p}(t)$ for simplicity.
Note that the BSplineManifold
objects are callable, and the arguments will be checked if it fits in the domain of BSplineSpace
.
BasicBSpline.BSplineManifold
— TypeConstruct 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 range.
[...]
If you need extension of BSplineManifold
or don't need the arguments check, you can call unbounded_mapping
.
BasicBSpline.unbounded_mapping
— Functionunbounded_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 range.
[...]
unbounded_mapping(M,t...)
is a little bit faster than M(t...)
because it does not check the domain.
B-spline curve
## 1-dim B-spline manifold
p = 2 # degree of polynomial
k = KnotVector(1:12) # knot vector
P = BSplineSpace{p}(k) # B-spline space
a = [SVector(i-5, 3*sin(i^2)) for i in 1:dim(P)] # control points
M = BSplineManifold(a, P) # Define B-spline manifold
plot(M)
B-spline surface
## 2-dim B-spline manifold
p = 2 # degree of polynomial
k = KnotVector(1:8) # knot vector
P = BSplineSpace{p}(k) # B-spline space
rand_a = [SVector(rand(), rand(), rand()) for i in 1:dim(P), j in 1:dim(P)]
a = [SVector(2*i-6.5, 2*j-6.5, 0) for i in 1:dim(P), j in 1:dim(P)] + rand_a # random generated control points
M = BSplineManifold(a,(P,P)) # Define B-spline manifold
plot(M)
Affine commutativity
Let $T$ be a affine transform $V \to W$, then the following equality holds.
\[T(\bm{p}(t; \bm{a})) =\bm{p}(t; T(\bm{a}))\]
Fixing arguments (currying)
Just like fixing first index such as A[3,:]::Array{1}
for matrix A::Array{2}
, fixing first argument M(4.3,:)
will create BSplineManifold{1}
for B-spline surface M::BSplineManifold{2}
.
julia> p = 2;
julia> k = KnotVector(1:8);
julia> P = BSplineSpace{p}(k);
julia> a = [SVector(5i+5j+rand(), 5i-5j+rand(), rand()) for i in 1:dim(P), j in 1:dim(P)];
julia> M = BSplineManifold(a,(P,P));
julia> M isa BSplineManifold{2}
true
julia> M(:,:) isa BSplineManifold{2}
true
julia> M(4.3,:) isa BSplineManifold{1} # Fix first argument
true
plot(M)
plot!(M(4.3,:), linewidth = 5, color=:cyan)
plot!(M(4.4,:), linewidth = 5, color=:red)
plot!(M(:,5.2), linewidth = 5, color=:green)