UnivariateSplines

Index

IgaBase.dimsplinespaceMethod
dimsplinespace(p, U)

Compute the dimension of the spline space defined by degree $p$ and knotsvector $U$.

source
IgaBase.refineMethod
refine(p::Degree, R::AbstractRefinement)

Change in polynomial degree based on the refinement type.

source
IgaBase.refineMethod
refine(p::KnotVector, R::AbstractRefinement)

Change in knot vector based on the refinement type.

source
UnivariateSplines.approximate_collocation_inverseMethod
approximate_collocation_inverse(p, U)
approximate_collocation_inverse(p, U, k1)

Compute an approximate inverse $\mathsf{A} ≈ \mathsf{B}^{-1}$, where $\mathsf{B}$ is the consistent B-spline interpolation matrix. The approximation is of order $k_1 \leq p+1$. The quasi-interpolant is designed to reproduce polynomials

\[ \sum_{i=1}^n \mu_i(x^k) B_{i,p}(x) = x^k  \quad \text{for } k=0,1,...,k_1-1\]

and is based on the following paper

T. Lyche and L. L. Schumaker. Local spline approximation methods. Journal of Approximation Theory, 15(4):294-325, Dec. 1975.

Examples:

julia> p  = Degree(3);

julia> U  = KnotVector([0.0,1.0,1.5,2.5,3.0], [p+1,1,1,2,p+1]);

julia> y  = grevillepoints(p, U);

The approximate collocation inverse, using full approximation order $k_1=p+1$ and the consistent collocation matrix are computed as follows

julia> A = approximate_collocation_inverse(p, U, p+1);

julia> B = bspline_interpolation_matrix(p, U, y, 1)[1];

It can be verified that the quasi-interpolant reproduces constants, linears, quadratic and cubic polynomials

julia> B * A * (y.^0) ≈ y.^0
true

julia> B * A * (y.^1) ≈ y.^1
true

julia> B * A * (y.^2) ≈ y.^2
true

julia> B * A * (y.^3) ≈ y.^3
true
source
UnivariateSplines.approximate_l2_inverseMethod
approximate_l2_inverse(p, U, k1)

Compute an approximate inverse $\mathsf{S} ≈ \mathsf{M}^{-1}$, where $\mathsf{M}$ is the consisten B-spline mass matrix. The approximation is of order $k_1 \leq p+1$. The $L^2$ quasi-interpolant is designed to reproduce polynomials

\[ \sum_{i=1}^n \mu_i(x^k) B_{i,p}(x) = x^k  \quad \text{for } k=0,1,...,k_1-1\]

and is based on the following paper

Chui, Charles K., Wenjie He, and Joachim Stöckler. "Nonstationary tight wavelet frames, I: Bounded intervals." Applied and Computational Harmonic Analysis 17.2 (2004): 141-197.

source
UnivariateSplines.bezier_extraction_operatorMethod
bezier_extraction_operator(p, U)

Compute the Bézier extraction operators corresponding to a B-spline basis of polynomial degree $p$ and knot vector $U$. The output is a 3-dimensional array $\mathsf{C} \in \mathbb{R}^{(p+1) \times (p+1) \times n_{el}}$.

The implementation is based on the following paper

Borden, Michael J., et al. "Isogeometric finite element data structures based on Bézier extraction of NURBS." International Journal for Numerical Methods in Engineering 87.1‐5 (2011): 15-47.

Examples:

julia> p = Degree(2); U = KnotVector([0.0,0.0,0.0,1.0,3.0,3.0,4.0,4.0,4.0]);

julia> bezier_extraction_operator(p, U)
3×3×3 Array{Float64, 3}:
[:, :, 1] =
 1.0  0.0       0.0
 0.0  1.0       0.0
 0.0  0.666667  0.333333

[:, :, 2] =
 0.666667  0.333333  0.0
 0.0       1.0       0.0
 0.0       0.0       1.0

[:, :, 3] =
 1.0  0.0  0.0
 0.0  1.0  0.0
 0.0  0.0  1.0
source
UnivariateSplines.bspline_interpolation_matrixFunction
bspline_interpolation_matrix(p, U, i, x, nout=1)

Given the polynomial degree $p$ and knot-vector $U$ compute the collocation matrices of the B-spline basis functions sampled at a set of points $x$. The output is a set of sparse matrices whose rows correspond to the evaluation points and whose columns are associated with the B-spline degrees of freedom.

Examples:

julia> p = Degree(2);

julia> U = KnotVector([0.0,1.0,1.5,2.5,3.0], [3,1,1,2,3]);

julia> x = IncreasingVector([0.5,1.5,2.75]);

julia> B = bspline_interpolation_matrix(p, U, x, p+1);

The interpolation matrix for quadratic B-splines is

julia> B[1]
3×7 SparseArrays.SparseMatrixCSC{Float64, Int64} with 9 stored entries:
 0.25  0.583333  0.166667   ⋅         ⋅     ⋅    ⋅ 
  ⋅     ⋅        0.666667  0.333333  0.0    ⋅    ⋅ 
  ⋅     ⋅         ⋅         ⋅        0.25  0.5  0.25

The interpolation matrix for linear and constant B-splines is given by

julia> B[2]
3×6 SparseArrays.SparseMatrixCSC{Float64, Int64} with 6 stored entries:
 0.5  0.5   ⋅    ⋅    ⋅    ⋅ 
  ⋅    ⋅   1.0  0.0   ⋅    ⋅ 
  ⋅    ⋅    ⋅    ⋅   0.5  0.5

julia> B[3]
3×5 SparseArrays.SparseMatrixCSC{Float64, Int64} with 3 stored entries:
 1.0   ⋅    ⋅    ⋅    ⋅ 
  ⋅    ⋅   1.0   ⋅    ⋅ 
  ⋅    ⋅    ⋅    ⋅   1.0

It can be easily verified that the B-splines satisfy a partition of unity

julia> sum(B[1],dims=2)
3×1 Matrix{Float64}:
 1.0
 1.0
 1.0

julia> sum(B[2],dims=2)
3×1 Matrix{Float64}:
 1.0
 1.0
 1.0

julia> sum(B[3],dims=2)
3×1 Matrix{Float64}:
 1.0
 1.0
 1.0

The collocation matrices can be used to solve an interpolation problem or a PDE using the collocation method. For example, the interpolation matrix evaluated at the Greville points can be computed as follows

julia> x = grevillepoints(p, U);

julia> B = bspline_interpolation_matrix(p, U, x, p+1)[1]
7×7 SparseArrays.SparseMatrixCSC{Float64, Int64} with 21 stored entries:
 1.0   0.0        0.0        ⋅          ⋅     ⋅    ⋅ 
 0.25  0.583333   0.166667   ⋅          ⋅     ⋅    ⋅ 
  ⋅    0.0833333  0.833333  0.0833333   ⋅     ⋅    ⋅ 
  ⋅     ⋅         0.166667  0.583333   0.25   ⋅    ⋅ 
  ⋅     ⋅          ⋅         ⋅         1.0   0.0  0.0
  ⋅     ⋅          ⋅         ⋅         0.25  0.5  0.25
  ⋅     ⋅          ⋅         ⋅         0.0   0.0  1.0

The matrix $\mathsf{B}$ can then be used to interpolate a univariate function. For example, up to quadratic polynomials are exactly reproduced

# sample grid to test polynomial reproduction
julia> y = global_insert(x, 4);

julia> C = bspline_interpolation_matrix(p, U, y, 1)[1];

julia> C * (B \ x.^0) ≈ y.^0
true

julia> C * (B \ x.^1) ≈ y.^1
true

julia> C * (B \ x.^2) ≈ y.^2
true

Here matrix $\mathsf{C}$ evaluates the B-spline at a refined partition $y$.

source
UnivariateSplines.bsplinebasisfunsFunction
bsplinebasisfuns(p, U, i::Integer, x::T, nout=1)
bsplinebasisfuns(p, U, x::T, nout=1)

Compute the $p+1$ non-zero B-spline basis-functions at site $x \in [U_{i},U_{i+1})$. The output is a matrix where column $j$ (rows $1:p+2-j$), for $j=1,...,n_{out}$ correspond to the $p+2-j$ non-zero B-spline basis functions of degree $p+1-j$. The remaining entries correspond to supports. If the active knot-span is not provided then it will be computed.

Examples:

julia> p = Degree(2);

julia> U = KnotVector([0.0,0.0,0.0,1.0,2.0,2.0,2.0]);

julia> x = 1.5
1.5

julia> span = findspan(p, U, x)
4

julia> B = bsplinebasisfuns(p, U, span, x, p+1)
3×3 Matrix{Float64}:
 0.125  0.5  1.0
 0.625  0.5  1.0
 0.25   1.0  2.0

The knot span index can be omited, in which case it is computed on the fly.

julia> bsplinebasisfuns(p, U, x, p+1)
3×3 Matrix{Float64}:
 0.125  0.5  1.0
 0.625  0.5  1.0
 0.25   1.0  2.0

By prescribing a vector of $m$ sites $x \in \mathbb{R}^m$ we can compute the B-splines at more than one point. The function can be called with or without prescribing the knot-span index

bsplinebasisfuns(p, U, span::Vector{Integer}, u::Vector{T}, nout=1)
bsplinebasisfuns(p, U, u::SortedSequence{T}, nout=1)

The output is a 3-dimensional array $B \in \mathbb{R}^{(p+1) \times m \times n_{out}}$.

Examples:

julia> x = IncreasingVector([0.5, 1.5]);

julia> span = findspan(p, U, x)
2-element NonDecreasingVector{Int64}:
 3
 4

julia> bsplinebasisfuns(p, U, span, x, p+1)
3×2×3 Array{Float64, 3}:
[:, :, 1] =
 0.25   0.125
 0.625  0.625
 0.125  0.25

[:, :, 2] =
 0.5  0.5
 0.5  0.5
 2.0  1.0

[:, :, 3] =
 1.0  1.0
 1.0  1.0
 1.0  2.0

Alternatively, the non-vanishing B-splines up to and including degree $p$ can be computed at $n_{el}$ non-zero knot-spans (elements) with $m$ sites each, that is, $x \in \mathbb{R}^{m \times n_{el}}$. Again, the function can be called with or without prescribing the knot-span index

bsplinebasisfuns(p, U, span::Vector{Int64}, u::Matrix{T}, nout=1)
bsplinebasisfuns(p, U, u::Matrix{T}, nout=1)

The output is a 4-dimensional array $B \in \mathbb{R}^{(p+1) \times m \times n_{el} \times n_{out}}$.

Examples:

julia> x = [0.25 1.25;
            0.75 1.75];

julia> bsplinebasisfuns(p, U, x, p+1)
3×2×2×3 Array{Float64, 4}:
[:, :, 1, 1] =
 0.5625   0.0625
 0.40625  0.65625
 0.03125  0.28125

[:, :, 2, 1] =
 0.28125  0.03125
 0.65625  0.40625
 0.0625   0.5625

[:, :, 1, 2] =
 0.75  0.25
 0.25  0.75
 2.0   2.0

[:, :, 2, 2] =
 0.75  0.25
 0.25  0.75
 1.0   1.0

[:, :, 1, 3] =
 1.0  1.0
 1.0  1.0
 1.0  1.0

[:, :, 2, 3] =
 1.0  1.0
 1.0  1.0
 2.0  2.0
source
UnivariateSplines.ders_bspline_interpolation_matrixFunction
ders_bspline_interpolation_matrix(p, U, i, x, nout=1)

Given the polynomial degree $p$ and knot-vector $U$ compute the collocation matrices of the first $n_{out}$ derivatives of the B-spline basis functions sampled at a set of points $x$. The output is a set of sparse matrices whose rows correspond to the evaluation points and whose columns are associated with the B-spline degrees of freedom.

Examples:

julia> p = Degree(2);

julia> U = KnotVector([0.0,1.0,1.5,2.5,3.0], [3,1,1,2,3]);

julia> x = IncreasingVector([0.5,1.5,2.75]);

julia> B = ders_bspline_interpolation_matrix(p, U, x, p+1);

The output is a set of SparseMatrixCSC matrices. For example the B-splines evaluate at points $x$ are

julia> B[1]
3×7 SparseArrays.SparseMatrixCSC{Float64, Int64} with 9 stored entries:                       
 0.25  0.583333  0.166667   ⋅         ⋅     ⋅    ⋅                                            
  ⋅     ⋅        0.666667  0.333333  0.0    ⋅    ⋅                                            
  ⋅     ⋅         ⋅         ⋅        0.25  0.5  0.25                                          

Derivatives of quadratic B-splines evaluated at $x$

julia> B[2]
3×7 SparseArrays.SparseMatrixCSC{Float64, Int64} with 9 stored entries:
 -1.0  0.333333   0.666667   ⋅         ⋅    ⋅    ⋅ 
   ⋅    ⋅        -1.33333   1.33333   0.0   ⋅    ⋅ 
   ⋅    ⋅          ⋅         ⋅       -2.0  0.0  2.0

Second derivatives of quadratic B-splines evaluated at $x$

julia> B[3]
3×7 SparseArrays.SparseMatrixCSC{Float64, Int64} with 9 stored entries:
 2.0  -3.33333  1.33333    ⋅        ⋅      ⋅    ⋅ 
  ⋅     ⋅       1.33333  -3.33333  2.0     ⋅    ⋅ 
  ⋅     ⋅        ⋅         ⋅       8.0  -16.0  8.0
source
UnivariateSplines.dersbsplinebasisfunsFunction
dersbsplinebasisfuns(p, U, span::Int, u::T, nout=1)
dersbsplinebasisfuns(p, U, u::T, nout=1)

Compute the non-vanishing B-spline basis-functions of degree $p$ and their $n_{out}-1$'th order derivatives at the site $x$. The output is a matrix $B \in \mathbb{R}^{(p+1) \times n_{out}}$ where the $i$'th row corresponds to the $i$'th non-vanishing B-spline at $x$ and the $j$'th column corresponds to the $j-1$'th derivative.

The implementation of this function is based on

Piegl, Les, and Wayne Tiller. The NURBS book. Springer Science & Business Media, 2012.

Examples:

julia> p = Degree(2);

julia> U = KnotVector([0.0,0.0,0.0,1.0,2.0,2.0,2.0]);

julia> x = 1.5;

julia> span = findspan(p, U, x)
4

julia> dersbsplinebasisfuns(p, U, span, x, p+1)
3×3 Matrix{Float64}:
 0.125  -0.5   1.0
 0.625  -0.5  -3.0
 0.25    1.0   2.0

dersbsplinebasisfuns can also be called without providing the span-indices, in which case they will be calculated on the fly.

julia> dersbsplinebasisfuns(p, U, x, p+1)
3×3 Matrix{Float64}:
 0.125  -0.5   1.0
 0.625  -0.5  -3.0
 0.25    1.0   2.0

By prescribing a vector of $m$ sites $x \in \mathbb{R}^m$ we can compute the B-splines and their derivatives at more than one point. The function can be called with or without prescribing the knot-span index

dersbsplinebasisfuns(p, U, span::Vector{Integer}, u::Vector{T}, nout=1)
dersbsplinebasisfuns(p, U, u::SortedSequence{T}, nout=1)

The output is a 3-dimensional array $B \in \mathbb{R}^{(p+1) \times m \times n_{out}}$.

Examples:

julia> x = IncreasingVector([0.5, 1.5]);

julia> span = findspan(p, U, x)
2-element NonDecreasingVector{Int64}:
 3
 4

julia> dersbsplinebasisfuns(p, U, span, x, p+1)
3×2×3 Array{Float64, 3}:
[:, :, 1] =
 0.25   0.125
 0.625  0.625
 0.125  0.25

[:, :, 2] =
 -1.0  -0.5
  0.5  -0.5
  0.5   1.0

[:, :, 3] =
  2.0   1.0
 -3.0  -3.0
  1.0   2.0

Alternatively, the non-vanishing derivatives of the B-splines can be computed at $n_{el}$ non-zero knot-spans (elements) with $m$ sites each, that is, $x \in \mathbb{R}^{m \times n_{el}}$. Again, the function can be called with or without prescribing the knot-span index

dersbsplinebasisfuns(p, U, span::Vector{Int64}, u::Matrix{T}, nout=1)
dersbsplinebasisfuns(p, U, u::Matrix{T}, nout=1)

The output is a 4-dimensional array $B \in \mathbb{R}^{(p+1) \times m \times n_{el} \times n_{out}}$.

Examples:

julia> x = [0.25 1.25;
            0.75 1.75];

julia> dersbsplinebasisfuns(p, U, x, p+1)
3×2×2×3 Array{Float64, 4}:
[:, :, 1, 1] =
 0.5625   0.0625
 0.40625  0.65625
 0.03125  0.28125

[:, :, 2, 1] =
 0.28125  0.03125
 0.65625  0.40625
 0.0625   0.5625

[:, :, 1, 2] =
 -1.5   -0.5
  1.25  -0.25
  0.25   0.75

[:, :, 2, 2] =
 -0.75  -0.25
  0.25  -1.25
  0.5    1.5

[:, :, 1, 3] =
  2.0   2.0
 -3.0  -3.0
  1.0   1.0

[:, :, 2, 3] =
  1.0   1.0
 -3.0  -3.0
  2.0   2.0
source
UnivariateSplines.distribute_pointsMethod
distribute_points(S::SplineSpace, V::SplineSpace; add_boundary_points::Bool=false)

Compute distribution of quadrature points in weigthed quadrature with a test-space $V$ and target space for quadrature $S$. The nodes are distribute in such a way that all conditions for exact quadrature are satisfied with a minimum number of points.

The distribution of quadrature points in weighted quadrature is based on the following paper.

Hiemstra, René R., et al. "Fast formation and assembly of finite element matrices with application to isogeometric linear elasticity." Computer Methods in Applied Mechanics and Engineering 355 (2019): 234-260.

Examples:

julia> S = SplineSpace(2, IncreasingVector([0.0,1.0,2.0,3.0]), [3,2,2,3]);

julia> V = SplineSpace(2, IncreasingVector([0.0,1.0,2.0,3.0]), [3,1,1,3]);

julia> UnivariateSplines.distribute_points(S, V)
8-element IncreasingVector{Float64}:
 0.16666666666666666
 0.5
 0.8333333333333334
 1.25
 1.75
 2.1666666666666665
 2.5
 2.8333333333333335
source
UnivariateSplines.extraction_operatorMethod
extraction_operator(p::Degree, U::KnotVector; [cperiodic, cleft, cright])

Compute an extraction operator that creates a subspace of the SplineSpace with for example peroidicity built-in.

source
UnivariateSplines.findspanMethod
findspan(p, U, u)

Given polynomial degree $p$ and knotvector $U$ determine the knot span index $i$ of a point $x$ such that $u \in [U_i, U_{i+1})$.

Examples:

julia> p = Degree(2)
2

julia> U = KnotVector([0.0,1.0,1.5,2.5,3.0],[3,1,1,2,3])
10-element NonDecreasingVector{Float64}:
 0.0
 0.0
 0.0
 1.0
 1.5
 2.5
 2.5
 3.0
 3.0
 3.0

julia> x = 2.0
2.0

julia> span = findspan(p, U, x)
5

julia> U[span] ≤ x < U[span+1]
true

It is also possible to compute the knot-spans of all values in an 'IncreasingVector{T}'

julia> x = IncreasingVector([0.25, 0.75, 1.0, 1.25, 2.5, 2.75]);

julia> span = findspan(p, U, x)
6-element NonDecreasingVector{Int64}:
 3
 3
 4
 4
 7
 7
source
UnivariateSplines.grevillepointsMethod
grevillepoints(p, U, i)
grevillepoints(p, U)

Compute the $i$'th Greville Absissa corresponding to the $i$'th B-spline defined by degree 'p' and knotvector 'U' of spline space.

Examples:

julia> p = Degree(2);

julia> U = KnotVector([0.0,1.0,2.5,3.0], [3,1,2,1])
7-element NonDecreasingVector{Float64}:
 0.0
 0.0
 0.0
 1.0
 2.5
 2.5
 3.0

julia> grevillepoints(p, U, 2)
0.5

julia> grevillepoints(p, U)
4-element IncreasingVector{Float64}:
 0.0
 0.5
 1.75
 2.5
source
UnivariateSplines.h_refinement_operator!Method
h_refinement_operator!(p::Degree, U::KnotVector{T}, u::AbstractVector{T}) where {T<:Real}
h_refinement_operator!(p::Degree, U::KnotVector{T}, u::T) where {T<:Real}
h_refinement_operator!(p::Degree, U::KnotVector{T}, span::Integer, u::T) where {T<:Real}

Insert one or more knots $u$ into knotvector $U$ and output the transformation operator from the coarse to the refined space. The knotvector is updated in-place. The output is represented as a SparseMatrixCSC matrix.

Examples:

julia> p = Degree(2)
2

julia> U = KnotVector([0.0,1.0,2.0,3.0,4], [3,1,1,2,3])
10-element NonDecreasingVector{Float64}:
 0.0
 0.0
 0.0
 1.0
 2.0
 3.0
 3.0
 4.0
 4.0
 4.0

julia> C = h_refinement_operator!(p, U, [0.5, 1.5, 2.5, 3.5]); # sparse matrix

julia> Matrix(C)
11×7 Matrix{Float64}:
 1.0  0.0   0.0   0.0   0.0  0.0  0.0
 0.5  0.5   0.0   0.0   0.0  0.0  0.0
 0.0  0.75  0.25  0.0   0.0  0.0  0.0
 0.0  0.25  0.75  0.0   0.0  0.0  0.0
 0.0  0.0   0.75  0.25  0.0  0.0  0.0
 0.0  0.0   0.25  0.75  0.0  0.0  0.0
 0.0  0.0   0.0   0.5   0.5  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   0.5  0.5  0.0
 0.0  0.0   0.0   0.0   0.0  0.5  0.5
 0.0  0.0   0.0   0.0   0.0  0.0  1.0

julia> @show U;
U = [0.0, 0.0, 0.0, 0.5, 1.0, 1.5, 2.0, 2.5, 3.0, 3.0, 3.5, 4.0, 4.0, 4.0]
source
UnivariateSplines.normalization_weightsMethod
normalization_weights(X)

Determine a set of normalization weights that result in a more uniform set oftype( weighted quadrature weights when solving for the least norm solution. )

source
UnivariateSplines.nquadpointsMethod
nquadpoints(S::SplineSpace, V::SplineSpace; add_boundary_points=false, add_additional_points=0, min_points=1)

Determine the minimum number of quadrature points corresponding to a weighted quadrature rule with test-space $V$ and target-space $S$.

Examples:

julia> S = SplineSpace(2, IncreasingVector([0.0,1.0,2.0,3.0]), [3,2,2,3]);

julia> V = SplineSpace(2, IncreasingVector([0.0,1.0,2.0,3.0]), [3,1,1,3]);

julia> A = UnivariateSplines.nquadpoints(S, V; add_boundary_points=false)
3-element Vector{Int64}:
 3
 2
 3
source
UnivariateSplines.onebasisfunevalMethod
onebasisfuneval(U, u)

Compute value of a single B-spline basis-function defined by local knotvector 'U' evaluated at point 'x'.

The implementation of this function is taken from

Piegl, Les, and Wayne Tiller. The NURBS book. Springer Science & Business Media, 2012.

source
UnivariateSplines.refinement_operatorMethod
refinement_operator(p::Degree, U::KnotVector, method::AbstractRefinement)

Compute the two_scale_operator given p and U and the applicable refinement strategy: IgaBase.hRefinement(), IgaBase.pRefinement(), IgaBase.kRefinement(). See IgaBase package for definition of refinement methods.

source
UnivariateSplines.table_required_pointsMethod
table_required_points(S, V)

Determine the number of exactness conditions in subsets of elements. The output is a matrix $\mathsf{A}$ where $\mathsf{A}_{ij}$ denotes the minimum number of points required in interval $i:j$ to perform exact weighted quadrature.

Examples:

julia> S = SplineSpace(2, IncreasingVector([0.0,1.0,2.0,3.0]), [3,2,2,3]);

julia> V = SplineSpace(2, IncreasingVector([0.0,1.0,2.0,3.0]), [3,1,1,3]);

julia> A = UnivariateSplines.table_required_points(S, V)
3×3×1 Array{Int64, 3}:
[:, :, 1] =
 3  5  7
 2  5  0
 3  0  0
source
UnivariateSplines.two_scale_operatorMethod
two_scale_operator(p::Degree, U::KnotVector{T}, q::Degree, V::KnotVector{T}) where {T<:Real}

The two_scale_operator encodes the main logic for all spline refinemement operations.

The arguments p and U correspond to the current spline space and q and V are the "target" degree knot vector.

source
UnivariateSplines.BsplineSupportType
BsplineSupport(S::SplineSpace)

Iterator that outputs the supporting elements of each B-spline basis-function as a UnitRange{Int64}.

Example:

julia> S = SplineSpace(2, [0.0,2.0,3.0], [3,1,3]);

julia> for α in BsplineSupport(S)
           @show α
       end
α = 1:1
α = 1:2
α = 1:2
α = 2:2
source
UnivariateSplines.GeneralizedGaussruleType
GeneralizedGaussrule(p, r, nel, a, b)

Compute a generalized Gaussian quadrature rule for a target spline-space $\mathbb{S}^p_{r}(a,b)$ defined in a uniform partition with $n_{el}$ elements.

These rules are made available in the following paper

Hiemstra, René R., et al. "Optimal and reduced quadrature rules for tensor product and hierarchically refined splines in isogeometric analysis." Computer Methods in Applied Mechanics and Engineering 316 (2017): 966-1004.

source
UnivariateSplines.SpanIndexType
SpanIndex(S::SplineSpace)

Iterator that outputs the span index corresponding to each element in the partition

Example:

julia> S = SplineSpace(2, [0.0,1.0,2.0,3.0,4.0], [3,1,2,1,3]);

julia> for s in SpanIndex(S)
           @show s
       end
s = 3
s = 4
s = 6
s = 7
source
UnivariateSplines.SplineSpaceType
SplineSpace(p, U)
SplineSpace(p, x, m)
SplineSpace(p, Interval(a, b), num_elements)

Definition of a spline space by means of the polynomial degree $p$, a sequence of break-points $x$ and the knot-multiplicity $m$ or by prescribing the knot-vector $U$.

Examples:

julia> S = SplineSpace(2, [0.0,2.0,3.0], [3,1,3])
SplineSpace(degree = 2, interval = [0.0, 3.0], dimension = 4)

julia> dimsplinespace(S)
4

julia> Degree(S)
2

julia> KnotVector(S)
7-element NonDecreasingVector{Float64}:
 0.0
 0.0
 0.0
 2.0
 3.0
 3.0
 3.0
source
UnivariateSplines.SupportType
BsplineSupport(S::SplineSpace)

Iterator that outputs the supporting elements of each B-spline basis-function as a UnitRange{Int64}.

source
UnivariateSplines.WeightedQuadruleType
WeightedQuadrule{T<:Real} <: AbstractQuadrule{1}

Weighted quadrature rules are test function specific quadrature rules that exactly integrate all functions in a targetspace. The rules are computed in high precision and are accurate up to 16 digits.

source