B-spline manifold

Setup

using BasicBSpline
using StaticArrays
using StaticArrays
using Plots

Multi-dimensional B-spline

Thm. Basis of tensor product of B-spline spaces

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)\]

The next plot shows $B_{(3,p^1,k^1)} \otimes B_{(4,p^2,k^2)}$ basis function.

# Define shape
k1 = KnotVector([0.0, 1.5, 2.5, 5.5, 8.0, 9.0, 9.5, 10.0])
k2 = knotvector"31 2  121"
P1 = BSplineSpace{3}(k1)
P2 = BSplineSpace{2}(k2)

# Choose indices
i1 = 3
i2 = 4

# Visualize basis functions
plotly()
plot(P1; plane=:xz, label="P1", color=:red)
plot!(P2; plane=:yz, label="P2", color=:green)
plot!(k1; plane=:xz, label="k1", color=:red, markersize=2)
plot!(k2; plane=:yz, label="k2", color=:green, markersize=2)
xs = range(0,10,length=100)
ys = range(1,9,length=100)
surface!(xs, ys, bsplinebasis.(P1,i1,xs') .* bsplinebasis.(P2,i2,ys))

Higher dimensional tensor products $\mathcal{P}[p^1,k^1]\otimes\cdots\otimes\mathcal{P}[p^d,k^d]$ are defined similarly.

Definition

B-spline manifold is a parametric representation of a shape.

Def. B-spline manifold

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 parametrization:

\[\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.

P = BSplineSpace{2}(KnotVector([0,0,0,1,1,1]))
a = [SVector(1,0), SVector(1,1), SVector(0,1)]  # `length(a) == dim(P)`
M = BSplineManifold(a, P)
M(0.4)  # Calculate `sum(a[i]*bsplinebasis(P, i, 0.4) for i in 1:dim(P))`
2-element StaticArraysCore.SVector{2, Float64} with indices SOneTo(2):
 0.84
 0.64

If you need extension of BSplineManifold or don't need the arguments check, you can call unbounded_mapping.

julia> M(0.4)2-element StaticArraysCore.SVector{2, Float64} with indices SOneTo(2):
 0.84
 0.64
julia> unbounded_mapping(M, 0.4)2-element StaticArraysCore.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.
julia> unbounded_mapping(M, 1.2)2-element StaticArraysCore.SVector{2, Float64} with indices SOneTo(2): -0.4399999999999999 0.9600000000000001

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)

Paraboloid

plotly()
p = 2
k = KnotVector([-1,-1,-1,1,1,1])
P = BSplineSpace{p}(k)
a = [SVector(i,j,2i^2+2j^2-2) for i in -1:1, j in -1:1]
M = BSplineManifold(a,P,P)
plot(M)

Hyperbolic paraboloid

plotly()
a = [SVector(i,j,2i^2-2j^2) for i in -1:1, j in -1:1]
M = BSplineManifold(a,P,P)
plot(M)

B-spline solid

k1 = k2 = KnotVector([0,0,1,1])
k3 = UniformKnotVector(-6:6)
P1 = BSplineSpace{1}(k1)
P2 = BSplineSpace{1}(k2)
P3 = BSplineSpace{3}(k3)
e₁(t) = SVector(cos(t),sin(t),0)
e₂(t) = SVector(-sin(t),cos(t),0)
a = cat([[e₁(t)*i+e₂(t)*j+SVector(0,0,t) for i in 0:1, j in 0:1] for t in -2:0.5:2]..., dims=3)
M = BSplineManifold(a,(P1,P2,P3))
plot(M; colorbar=false)

Affine commutativity

Thm. 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 argumenttrue
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)