API Reference

FiniteElementHermite

FiniteElementHermite.D_νeFunction
D_νe(x::Real, k::Integer, ρ::AbstractVector{<:Real})

Returns the derivative of the even basis element centered at ρ[k] at x

source
FiniteElementHermite.D_νoFunction
D_νo(x::Real, k::Integer, ρ::AbstractVector{<:Real})

Returns the derivative of the odd basis element centered at ρ[k] at x

source
FiniteElementHermite.FEFunction
FE(x::AbstractVector{<:Real}, y::AbstractVector{<:Real})

Returns finite element representation of data (x, y), fitting local quadratics to determine derivative

source
FE(grid::AbstractVector{<:Real}, xy::Tuple{<:AbstractVector{<:Real},<:AbstractVector{<:Real}})

Returns finite element representation of data (x, y), mapped to grid

source
FiniteElementHermite.FE_repType
FE_rep(x::S, coeffs::T) where {S<:AbstractVector{<:Real},T<:AbstractVector{<:Real}}

Create Hermite-cubic finite-element representation on grid x with cofficients coeffs coeffs[2k] is the value of the function at x[k] coeffs[2k-1] is the derivative of the function at x[k]

An FE_rep is callable, giving the value of the finite-element representation, for example:

Y = FE_rep(x, coeffs)
Y(a) # gives value of the finite-element representation at x=a
source
(Y::FE_rep)(x::Real)

Functor for FE_rep, giving the value of the finite-element representation at x

source
FiniteElementHermite.I_νeFunction
I_νe(x::Real, k::Integer, ρ::AbstractVector{<:Real})

Returns the integral of the even basis element centered at ρ[k] from 0 to x

source
FiniteElementHermite.I_νoFunction
I_νo(x::Real, k::Integer, ρ::AbstractVector{<:Real})

Returns the integral of the odd basis element centered at ρ[k] from 0 to x

source
FiniteElementHermite.compute_D_basesFunction
compute_D_bases(X::AbstractVector{<:Real}, x::Real)

For grid X, find the derivative of the four basis functions at x This can be used with evaluate or evaluate_inbounds to compute the derivative of multiple FE_reps efficiently if they share the same grid

Returns (k, D_nu_ou, D_nu_eu, D_nu_ol, D_nu_el) where D_nu_ou is derivative of odd basis element centered at X[k] D_nu_eu is derivative of even basis element centered at X[k] D_nu_ol is derivative of odd basis element centered at X[k+1] D_nu_el is derivative of even basis element centered at X[k+1]

source
FiniteElementHermite.compute_basesFunction
compute_bases(X::AbstractVector{<:Real}, x::Real)

For grid X, find the value of the four basis functions at x This can be used with evaluate or evaluate_inbounds to compute the value of multiple FE_reps efficiently if they share the same grid

Returns (k, nu_ou, nu_eu, nu_ol, nu_el) where nu_ou is odd basis element centered at X[k] nu_eu is even basis element centered at X[k] nu_ol is odd basis element centered at X[k+1] nu_el is even basis element centered at X[k+1]

source
FiniteElementHermite.compute_both_basesFunction
compute_both_bases(X::AbstractVector{<:Real}, x::Real)

For grid X, find the value and derivative of the four basis functions at x This can be used with evaluate or evaluate_inbounds to compute the value and derivative of multiple FE_rep efficiently if they share the same grid

Returns (k, nu_ou, nu_eu, nu_ol, nu_el, D_nu_ou, D_nu_eu, D_nu_ol, D_nu_el) where nu_ou is odd basis element centered at X[k] nu_eu is even basis element centered at X[k] nu_ol is odd basis element centered at X[k+1] nu_el is even basis element centered at X[k+1] D_nu_ou is derivative of odd basis element centered at X[k] D_nu_eu is derivative of even basis element centered at X[k] D_nu_ol is derivative of odd basis element centered at X[k+1] D_nu_el is derivative of even basis element centered at X[k+1]

source
FiniteElementHermite.dual_inner_productFunction
dual_inner_product(integrands::Function, k1, k2, ρ::AbstractVector{<:Real}, order::Union{Nothing,Integer}=nothing)

Compute inner product for two functions simultaneously (integrands returns two values) in region where the k1-th and k2-th finite elements of grid ρ overlap Uses Gauss-Legendre quadature of requested order if given

source
FiniteElementHermite.evaluateFunction
evaluate(Y::FE_rep, k::Integer, nu_ou::Real, nu_eu::Real, nu_ol::Real, nu_el::Real)

Evaluate FE_rep Y assuming basis function values/derivatives have been computed by compute_bases, compute_D_bases, or compute_both_bases

`nu_ou` is value or derivative of odd  basis element centered at `X[k]`
`nu_eu` is value or derivative of even basis element centered at `X[k]`
`nu_ol` is value or derivative of odd  basis element centered at `X[k+1]`
`nu_el` is value or derivative of even basis element centered at `X[k+1]`
source
FiniteElementHermite.evaluate_inboundsFunction
evaluate_inbounds(Y::FE_rep, k::Integer, nu_ou::T, nu_eu::T, nu_ol::T, nu_el::T) where {T<:Real}

Evaluate FE_rep Y without bounds checking, assuming basis function values/derivatives have been computed by compute_bases, compute_D_bases, or compute_both_bases

`nu_ou` is value or derivative of odd  basis element centered at `X[k]`
`nu_eu` is value or derivative of even basis element centered at `X[k]`
`nu_ol` is value or derivative of odd  basis element centered at `X[k+1]`
`nu_el` is value or derivative of even basis element centered at `X[k+1]`
source
FiniteElementHermite.get_quadratureFunction
get_quadrature(x::AbstractVector, order::Integer)

For finite-element grid x, give Gauss-Legendre quadrature integration of requested order Returns (xq, wq), quadrature points and quadrature weights

source
FiniteElementHermite.inner_productFunction
inner_product(nu1, k1::Integer, nu2, k2::Integer, ρ::AbstractVector{<:Real}, order::Union{Nothing,Integer}=nothing)

Compute inner product ∫dx nu1(x, k1, ρ) * nu2(x, k2, ρ) where nu1 and nu2 are basis functions Uses Gauss-Legendre quadature of requested order if given

source
inner_product(f, nu1, k1::Integer, nu2, k2::Integer, ρ::AbstractVector{<:Real}, order::Union{Nothing,Integer}=nothing)

Compute inner product ∫dx f(x) * nu1(x, k1, ρ) * nu2(x, k2, ρ) where nu1 and nu2 are basis functions Uses Gauss-Legendre quadature of requested order if given

source
inner_product(f, nu1, k1::Integer, nu2, k2::Integer, ρ::AbstractVector{<:Real}, order::Union{Nothing,Integer}=nothing)

Compute inner product ∫dx nu1(x, k1, ρ) * (f(x) * fnu2(x, k2, ρ) + g(x) * gnu2(x, k2, ρ)) where nu1, fnu2, and gnu2 are basis functions Uses Gauss-Legendre quadature of requested order if given

source
inner_product(f::Function, nu::Function, k::Integer, ρ::AbstractVector{<:Real}, order::Union{Nothing,Integer}=nothing)

Compute inner product ∫dx f(x) * nu(x, k, ρ) where nu is a basis function Uses Gauss-Legendre quadature of requested order if given

source
inner_product(integrand::Function, k1, k2, ρ::AbstractVector{<:Real}, order::Union{Nothing,Integer}=nothing)

Compute inner product for generic integrand for in region where the k1-th and k2-th finite elements of grid ρ overlap Uses Gauss-Legendre quadature of requested order if given

source
FiniteElementHermite.mass_matrixFunction
mass_matrix(N::Integer, ρ::AbstractVector{<:Real})

NxN matrix of inner products of all Hermite cubic finite elements on grid ρ Returns a BandedMatrix as only neighboring finite elements overlap N.B. - N should be determined from ρ and will be removed later

source
FiniteElementHermite.νeFunction
νe(x::Real, k::Integer, ρ::AbstractVector{<:Real})

Returns the value of the even basis element centered at ρ[k] at x

source
FiniteElementHermite.νoFunction
νo(x::Real, k::Integer, ρ::AbstractVector{<:Real})

Returns the value of the odd basis element centered at ρ[k] at x

source