BezierBernsteinMethods
This package implements Bezier-Bernstein methods on simplicial meshes.
References
- Kirby, R.C. Fast simplicial finite element algorithms using Bernstein polynomials. Numer. Math. 117, 631–652 (2011). https://doi.org/10.1007/s00211-010-0327-2
- Mark Ainsworth, Gaelle Andriamaro, and Oleg Davydov. Bernstein-Bézier Finite Elements of Arbitrary Order and Optimal Assembly Procedures. SIAM J. Sci. Comput. 33, 6 (November 2011), 3087–3109. https://doi.org/10.1137/11082539X
- Kirby, R.C. Low-Complexity Finite Element Algorithms for the de Rham Complex on Simplices. SIAM Journal on Scientific Computing 36:2, A846-A868 (2014). https://doi.org/10.1137/130927693
Index
BezierBernsteinMethods.BezierSimplex
BezierBernsteinMethods.LinearIndex
BezierBernsteinMethods.LinearIndex
BezierBernsteinMethods.MultiIndices
BezierBernsteinMethods.Simplex
BezierBernsteinMethods.Binomial
BezierBernsteinMethods.b2c
BezierBernsteinMethods.barycentriccoords
BezierBernsteinMethods.bernstein
BezierBernsteinMethods.bernsteinfuns
BezierBernsteinMethods.cartesiancoords
BezierBernsteinMethods.decasteljau!
BezierBernsteinMethods.degree_elevation_operator
BezierBernsteinMethods.grevillepoint
BezierBernsteinMethods.jacobian
BezierBernsteinMethods.kformbasisconv
BezierBernsteinMethods.linindex
BezierBernsteinMethods.mass_matrix
BezierBernsteinMethods.multiindex
BezierBernsteinMethods.step_decasteljau!
BezierBernsteinMethods.subsets
BezierBernsteinMethods.vol
IgaBase.dimension
IgaBase.dimsplinespace
BezierBernsteinMethods.Binomial
— MethodBinomial(Val(n), Val(k)) = binomial(n,k)
@generated
function for optimized computation of binomials.
BezierBernsteinMethods.b2c
— Methodb2c(simplex)
Return the mapping
(matrix) that maps barycentric coordinates to Cartesian coordinates.
BezierBernsteinMethods.barycentriccoords
— Methodbarycentriccoords(S::AbstractSimplex, x)
Compute the barycentric coordinates of a point, or matrix of points x
.
BezierBernsteinMethods.bernstein
— Methodbernstein(p::Degree, control::AbstractVector{T}, λ::BPoint{D,T})
Compute Bezier polynomial at barycentric point λ
with the degrees of freedom in control
.
BezierBernsteinMethods.bernsteinfuns
— Methodbernstein(X::MultiIndex{D}, λ::BPoint{D,T})
Compute single Bernstein basis function corresponding to multi-index X
at barycentric point λ
.
BezierBernsteinMethods.cartesiancoords
— Methodcartesiancoords(S::AbstractSimplex, λ)
Compute the Cartesian coordinates of a barycentric point or vector of barycentric points λ
.
BezierBernsteinMethods.decasteljau!
— Methoddecasteljau!(control::Vector{T}, λ::BPoint{D,T}, MultiIndices)
Classic DeCasteljau algorithm (can compute derivatives too!) on a (D-1)
-dimensional simplex.
BezierBernsteinMethods.degree_elevation_operator
— Functiondegree_elevation_operator(p::Degree, d::Dimension)
Return a matrix that maps degrees of freedom of a Bezier polynomial of degree p
to the degrees of freedom corresponding to its order elevated Bezier polynomial.
BezierBernsteinMethods.grevillepoint
— Methodgrevillepoint(X::MultiIndex)
Compute Greville point associated with multi-index label X
.
BezierBernsteinMethods.jacobian
— Methodjacobian(simplex)
Return the Jacobian mapping of a simplex.
BezierBernsteinMethods.kformbasisconv
— Methodkformbasisconv(n::Dimension, k::Form, r::Degree)
Compute the linear mapping from Bernstein polynomials to the space of k
-forms Pᵣ⁻Λᵏ(T)
. This linear mapping is encoded in Pattern{T}
. We refer to the paper. For more information, check out the periodic table of finite elements
BezierBernsteinMethods.linindex
— Methodlinindex(i, j, k...)
Return the linear index of a multi-index label.
BezierBernsteinMethods.mass_matrix
— Methodmass_matrix(S::AbstractSimplex{D}, p::Degree, q::Degree, T::Type{<:Real} = Type{Float64})
Compute Bernstein mass matrix for degree (p,q)
and dimension d
.
BezierBernsteinMethods.multiindex
— MethodMultiIndex(::LinearIndex{P<:Degree, D<:Dimension, I<:Integer})
Return the multi-index associated with linear index I
and simplicial domain of degree P
and dimension D
.
BezierBernsteinMethods.step_decasteljau!
— Methodstep_decasteljau!(control::Vector{T}, λ::BPoint{D,T}, MultiIndices)
Computes a single step of the recursive DeCasteljau algorithm on a (D-1)
-dimensional simplex.
BezierBernsteinMethods.subsets
— Methodsubsets(v::NTuple{n,Int}, k::Int)
Generate the binomial(n,k)
unique subsets of length k
from set v
and its complement.
Examples:
julia> α, β = subsets((1,2,3), 1);
julia> α
3-element Array{Tuple{Int64},1}:
(1,)
(2,)
(3,)
julia> β
3-element Array{Tuple{Int64,Int64},1}:
(2, 3)
(1, 3)
(1, 2)
BezierBernsteinMethods.vol
— Methodvol(simplex)
Compute the volume of an affine simplex.
IgaBase.dimension
— Methoddimension(simplex)
Return the dimension of a simplex.
IgaBase.dimsplinespace
— Methoddimsplinespace(p::Degree, d::Dimension)
Compute the dimension of the Bezier spline space on a simplex of degree p
and numer of vertices d
.
BezierBernsteinMethods.BezierSimplex
— Type BezierSimplex{D,T<:Real} <: AbstractSimplex{D}
Simplex of dimension D
with Vertices in real number type T
.
BezierBernsteinMethods.LinearIndex
— TypeLinearIndex(P::Degree, D::Dimension, I::Integer)
A linear index I
of a simplicial discretization with Degree P
and dimension D
.
BezierBernsteinMethods.LinearIndex
— MethodLinearIndex(::MultiIndex{X})
Compute the LinearIndex
associated with a MultiIndex
.
BezierBernsteinMethods.MultiIndices
— TypeMultiIndices(P::Degree, D::Dimension)
MultiIndices(D::Dimension)
MultiIndices(s::Simplex)
MultiIndices(s::BezierSimplex)
Returns an iterator corresponding to a Bezier simplex of polynomial degree P
and dimension D
.
BezierBernsteinMethods.Simplex
— Type Simplex{D,T<:Real} <: AbstractSimplex{D}
Simplex of dimension D
with Vertices in real number type T
.