Base.Matrix
— MethodMatrix(L::DiagonalOperator[, x])
Compute the matrix representation of the DiagonalOperator
L
, optionally weighted by the function whose expansion coefficients are given by x
.
Note that this is allocating and potentially inefficient (especially in the case of non-orthogonal bases), and should thus not be used in compute-intense code. Furthermore, if you have a Julia function f
corresponding to the diagonal operator, it is most likely more reasonable to take the R'*QuasiDiagonal(f.(x))*R
approach instead.
CompactBases.AbstractKnotSet
— TypeAbstractKnotSet{k,ml,mr,T}
Abstract base for B-spline knot sets. T
is the eltype
of the knot set, k
is the order of the piecewise polynomials (order = degree + 1
) and ml
and mr
are the knot multiplicities of the left and right endpoints, respectively.
CompactBases.ArbitraryKnotSet
— TypeArbitraryKnotset(k, t[, ml=k, mr=k])
Construct an order-k
arbitrary knot set, with the locations of the knots given by non-decreasing vector t
.
CompactBases.ArbitraryKnotSet
— TypeArbitraryKnotSet{k,ml,mr}(t)
An arbitrary knot set of order k
and left and right multiplicities of ml
and mr
, respectively. The knot set is specified by the non-decreasing vector t
; the same knot may appear multiple times, which influences the continuity of the B-splines at that location.
CompactBases.BSpline
— TypeBSpline(t, x, w, B, S)
Basis structure for the B-splines generated by the knot set t
. x
and w
are the associated quadrature roots and weights, respectively, the columns of B
correspond to the B-splines resolved on the quadrature roots, and S
is the (banded) B-spline overlap matrix.
CompactBases.BSpline
— MethodBSpline(t[, N])
Create the B-spline basis corresponding to the knot set t
. N
is the amount of Gauß–Legendre quadrature points per interval.
CompactBases.BSpline
— MethodBSpline(t[; k′=3])
Create the B-spline basis corresponding to the knot set t
. k′
is the highest polynomial order of operators for which it should be possible to compute the matrix elements exactly (via Gauß–Legendre quadrature). The default k′=3
corresponds to operators O(x²).
CompactBases.BasisOrthogonality
— TypeBasisOrthogonality
Base type for the orthogonality trait of a basis.
CompactBases.CoulombDerivative
— TypeCoulombDerivative(D, Z, ℓ)
Helper operator wrapping the Derivative
operator in the case of a Coulomb potential of charge Z
and centrifugal potential $ℓ(ℓ+1)/2r²$. Its main use is to correctly apply boundary conditions at $r=0$, when materializing derivatives.
CompactBases.Density
— TypeDensity
Type-alias for FunctionProduct
where the first function is conjugated, as is necessary in complex linear algebra, when computing mutual densities.
CompactBases.DiagonalOperator
— TypeDiagonalOperator
Helper structure that when acting on a function $f(x)$ produces the product $f(x)g(x)$, where $g(x)$ is e.g. a potential. This is similar to R'*QuasiDiagonal(g.(x))*R
in intent, but simplifies the case when $g(x)$ needs to be updated, without computing a lot overlap integrals. This is accomplished by computing the FunctionProduct
of the two functions directly.
CompactBases.DiagonalOperator
— MethodDiagonalOperator(f)
Construct DiagonalOperator
from a function expansion f
.
CompactBases.ExpKnotSet
— TypeExpKnotSet{k,ml,mr}(exponents, base, t, include0)
A knot set of order k
and left and right multiplicities of ml
and mr
, respectively, whose knots are exponentially distributed according to t = base .^ exponents
, optionally including 0 as the left endpoint.
CompactBases.ExpKnotSet
— MethodExpKnotSet(k, a, b, N[, ml=k, mr=k, base=10, include0=true])
Construct an order-k
knot spanning from base^a
to base^b
in N
intervals, optionally including 0 as the left endpoint.
Examples
julia> ExpKnotSet(2, -4, 2, 7)
10-element ExpKnotSet{2,2,2,Float64,StepRangeLen{Float64,Base.TwicePrecision{Float64},Base.TwicePrecision{Float64}},Array{Float64,1}}:
0.0
0.0
0.0001
0.001
0.01
0.1
1.0
10.0
100.0
100.0
CompactBases.FunctionProduct
— TypeFunctionProduct{Conjugated}
Helper object to compute the expansion coefficients of $ρ(x) \defd f^\circ(x)g(x)$, where $f^\circ$ denotes that $f$ may be conjugated, if so desired.
CompactBases.FunctionProduct
— MethodFunctionProduct{Conjugated}(f, g) where Conjugated
Construct a FunctionProduct
for computing the product of the two functions f
and g
.
CompactBases.FunctionProduct
— MethodFunctionProduct{Conjugated}(L, R[, T; w=one]) where {Conjugated,T}
Construct a FunctionProduct
for computing the product of two functions expanded over L
and R
, respectively:
\[h(x) = f^\circ(x)g(x)w(x)\]
where $w(x)$ is an optional weight function, over the basis of $g(x)$, via Vandermonde interpolation.
CompactBases.LinearKnotSet
— TypeLinearKnotSet(k, a, b, N[, ml=k, mr=k])
Construct an order-k
linear knot set spanning from a
to b
, with N
intervals.
Examples
julia> LinearKnotSet(2, 0, 1, 3)
6-element LinearKnotSet{2,2,2,Float64,StepRangeLen{Float64,Base.TwicePrecision{Float64},Base.TwicePrecision{Float64}}}:
0.0
0.0
0.3333333333333333
0.6666666666666666
1.0
1.0
CompactBases.LinearKnotSet
— TypeLinearKnotSet{k,ml,mr}(t)
A knot set of order k
and left and right multiplicities of ml
and mr
, respectively, whose knots are uniformly distributed according to the range t
. The interior basis functions are thus Cᵏ⁻²-continuous.
CompactBases.LinearOperator
— TypeLinearOperator(A, B⁻¹, temp)
A helper object used to apply the action of a linear operator, whose matrix representation is given by A
, in the basis whose metric inverse is B⁻¹
. For orthogonal bases (such as e.g. FiniteDifferences
and FEDVR
), only multiplication by A
is necessary, but non-orthonal bases (such as BSpline
) necessitates the application of the metric inverse as well.
CompactBases.LinearOperator
— MethodLinearOperator(A, R)
Construct the linear operator whose matrix representation in the basis R
is A
, automatically deducing if application of the metric inverse is also necessary.
CompactBases.NonOrthogonal
— TypeNonOrthogonal <: BasisOrthogonality
Concrete type for the trait of non-orthogonal bases, i.e. with non-diagonal metric.
CompactBases.Ortho
— TypeOrtho <: BasisOrthogonality
Intermediate type for the trait for all orthogonal bases.
CompactBases.Orthogonal
— TypeOrthogonal <: Ortho
Concrete type for the trait of orthogonal bases, i.e. with diagonal metric.
CompactBases.Orthonormal
— TypeOrthonormal <: Ortho
Concrete type for the trait of orthonormal bases, i.e. with metric I
.
CompactBases.ShiftAndInvert
— TypeShiftAndInvert(A, ::UniformScaling[, σ=0])
Construct the shifted-and-inverted operator corresponding to the eigenproblem $(\mat{A}-σ\mat{I})\vec{x} = \lambda\vec{x}$.
CompactBases.ShiftAndInvert
— TypeShiftAndInvert(A, R[, σ=0])
Construct the shifted-and-inverted operator corresponding to the eigenproblem $(\mat{A}-σ\mat{B})\vec{x} = \lambda\mat{B}\vec{x}$, where $\mat{B}$ is the suitable operator metric, depending on the AbstractQuasiMatrix
R
employed.
CompactBases.ShiftAndInvert
— TypeShiftAndInvert(A, B[, σ=0])
Construct the shifted-and-inverted operator corresponding to the eigenproblem $(\mat{A}-σ\mat{B})\vec{x} = \lambda\mat{B}\vec{x}$.
CompactBases.ShiftAndInvert
— TypeShiftAndInvert(A⁻¹, B, temp)
Represents a shifted-and-inverted operator, suitable for Krylov iterations when targetting an interior point of the eigenspectrum. A⁻¹
is a factorization of the shifted operator and B
is the metric matrix.
CompactBases.StaggeredFiniteDifferences
— TypeStaggeredFiniteDifferences
Staggered finite differences with variationally derived three-points stencils for the case where there is Dirichlet0 boundary condition at r = 0. Supports non-uniform grids, c.f.
- Krause, J. L., & Schafer, K. J. (1999). Control of THz Emission from Stark Wave Packets. The Journal of Physical Chemistry A, 103(49), 10118–10125. http://dx.doi.org/10.1021/jp992144
CompactBases.UnitVector
— TypeUnitVector{T}(N, k)
Helper vector type of length N
where the k
th element is one(T)
and all the others zero(T)
.
Base.copyto!
— Methodcopyto!(o::DiagonalOperator, diag::AbstractVector)
Update the DiagonalOperator
to represent multiplication by the function whose expansion coefficients are diag
.
Base.copyto!
— MethodBase.copyto!(ρ::FunctionProduct, f, g)
Update the FunctionProduct
ρ
from the functions f
and g
.
Base.copyto!
— Methodcopyto!(ρ::FunctionProduct{Conjugated}, f::AbstractVector, g::AbstractVector) where Conjugated
Update the FunctionProduct
ρ
from the vectors of expansion coefficients, f
and g
.
Base.first
— Methodfirst(t)
Return the first knot of t
.
Base.getindex
— Methodgetindex(t, i)
Return the i
th knot of the knot set t
, accounting for the multiplicities of the endpoints.
Examples
julia> LinearKnotSet(3, 0, 1, 3)[2]
0.0
julia> LinearKnotSet(3, 0, 1, 3, 1, 1)[2]
0.3333333333333333
Base.last
— Methodlast(t)
Return the last knot of t
.
Base.length
— Methodlength(t)
Return the number of knots of t
.
Examples
julia> length(LinearKnotSet(3, 0, 1, 3))
8
julia> length(LinearKnotSet(3, 0, 1, 3, 1, 1))
4
CompactBases.assert_multiplicities
— Methodassert_multiplicities(k,ml,mr,t)
Assert that the multiplicities at the endpoints, ml
and mr
, respectively, are consistent with the order k
. Also check that the amount of knots in the knot set t
are enough to support the requested order k
.
CompactBases.basis_overlap
— Methodbasis_overlap(A::FiniteDifferencesOrRestricted, B)
Transforming to finite-differences is simple: just evaluate the basis functions of the source basis on the grid points, possibly weighting them.
CompactBases.basis_overlap
— Methodbasis_overlap(A::BSplineOrRestricted, B::Union{BSplineOrRestricted,FEDVROrRestricted})
From a piecewise polynomial basis, i.e. B-splines or FEDVR, we can use Gaussian quadrature to compute the overlap matrix elements exactly. It is possible that we could derive better results, if we used the fact that FEDVR basis functions are orthogonal in the sense of the Gauss–Lobatto quadrature, but we leave that for later.
CompactBases.basis_overlap
— Methodbasis_overlap(A::FEDVROrRestricted, B)
Transforming to FEDVR is much the same as for finite-differences; evaluate at the quadrature nodes and scale by the weights. However, this has to be done per-element, so we feed it through the interpolation routines already setup for this.
CompactBases.basis_transform
— Methodbasis_transform(A, B)
Compute the (possibly dense) basis transform matrix to go from B
to A
, via basis_overlap
, dispatching on the orthogonality
trait of A
.
CompactBases.basis_transform
— Methodbasis_transform(A, B)
Compute the (possibly dense) basis transform matrix to go from B
to A
, via basis_overlap
for non-orthogonal basis A
:
\[(B^HB)^{-1}B^HA\]
CompactBases.basis_transform
— Methodbasis_transform(A, B)
Compute the (possibly dense) basis transform matrix to go from B
to A
, via basis_overlap
for orthogonal basis A
:
\[B^HA\]
CompactBases.basis_transform
— Methodbasis_transform(A::BSplineOrRestricted, B::AbstractFiniteDifferences)
This works via Vandermonde interpolation, which is potentially unstable. In this case, we do not provide basis_overlap
.
CompactBases.centers
— Methodcenters(B)
Return the locations of the mass centers of all basis functions of B
; for orthogonal bases such as finite-differences and FE-DVR, this is simply locs
, i.e. the location of the quadrature nodes.
CompactBases.change_interval!
— Methodchange_interval!(xs, ws, x, w[, a=0, b=1, γ=1])
Transform the Gaußian quadrature roots x
and weights w
on the elementary interval [-1,1]
to the interval [γ*a,γ*b]
and store the result in xs
and ws
, respectively. γ
is an optional root of unity, used to complex-rotate the roots (but not the weights).
CompactBases.deBoor
— MethoddeBoor(t, c, x[, i[, m=0]])
Evaluate the spline given by the knot set t
and the set of control points c
at x
using de Boor's algorithm. i
is the index of the knot interval containing x
. If m≠0
, calculate the m
th derivative at x
instead.
CompactBases.find_interval
— Methodfind_interval(t, x[, i=ml])
Find the interval in the knot set t
that includes x
, starting from interval i
(which by default is the first non-zero interval of the knot set). The search complexity is linear, but by storing the result and using it as starting point for the next call to find_interval
, the knot set need only be traversed once.
CompactBases.findelement
— Functionfindelement(B::FEDVR, k[, i=1, m=k])
Find the finite-element of B
that contains the k
th basis function, optionally starting the search from element i
and basis function m
of that element.
CompactBases.lagrangeder!
— Methodlagrangeder!(xⁱ, m, L′)
Calculate the derivative of the Lagrange interpolating polynomial Lⁱₘ(x), given the roots xⁱ
, at the roots, and storing the result in L′
.
∂ₓ Lⁱₘ(xⁱₘ,) = (xⁱₘ-xⁱₘ,)⁻¹ ∏(k≠m,m′) (xⁱₘ,-xⁱₖ)/(xⁱₘ-xⁱₖ), m≠m′, [δ(m,n) - δ(m,1)]/2wⁱₘ, m=m′
Eq. (20) Rescigno2000
CompactBases.leftmultiplicity
— Methodleftmultiplicity(t)
Return the multiplicity of the left endpoint.
CompactBases.lgwt
— Methodlgwt(t, N) -> (x,w)
Generate the N
Gauß–Legendre quadrature roots x
and associated weights w
, with respect to the B-spline basis generated by the knot set t
.
Examples
julia> CompactBases.lgwt(LinearKnotSet(2, 0, 1, 3), 2)
([0.0704416, 0.262892, 0.403775, 0.596225, 0.737108, 0.929558], [0.166667, 0.166667, 0.166667, 0.166667, 0.166667, 0.166667])
julia> CompactBases.lgwt(ExpKnotSet(2, -4, 2, 7), 2)
([2.11325e-5, 7.88675e-5, 0.000290192, 0.000809808, 0.00290192, 0.00809808, 0.0290192, 0.0809808, 0.290192, 0.809808, 2.90192, 8.09808, 29.0192, 80.9808], [5.0e-5, 5.0e-5, 0.00045, 0.00045, 0.0045, 0.0045, 0.045, 0.045, 0.45, 0.45, 4.5, 4.5, 45.0, 45.0])
CompactBases.local_step
— Methodlocal_step(B, i)
The step size around grid point i
. For uniform grids, this is equivalent to step
.
CompactBases.metric
— Methodmetric(B)
Returns the metric or mass matrix of the basis B
, equivalent to S=B'B
.
CompactBases.metric_shape
— Methodmetric_shape(B)
Returns the shape of the metric or mass matrix of the basis B
, i.e. Diagonal
, BandedMatrix
, etc.
CompactBases.nonempty_intervals
— Methodnonempty_intervals(t)
Return the indices of all intervals of the knot set t
that are non-empty.
Examples
julia> nonempty_intervals(ArbitraryKnotSet(3, [0.0, 1, 1, 3, 4, 6], 1, 3))
4-element Array{Int64,1}:
1
3
4
5
CompactBases.num_quadrature_points
— Methodnum_quadrature_points(k, k′)
The number of quadrature points needed to exactly compute the matrix elements of an operator of polynomial order k′
with respect to a basis of order k
.
CompactBases.numfunctions
— Methodnumfunctions(t)
Returns the number of basis functions generated by knot set t
.
Examples
julia> numfunctions(LinearKnotSet(3, 0, 1, 2))
4
julia> numfunctions(LinearKnotSet(5, 0, 1, 2))
6
CompactBases.numintervals
— Methodnuminterval(t)
Returns the number of intervals generated by the knot set t
.
Examples
julia> numintervals(LinearKnotSet(3, 0, 1, 2))
2
CompactBases.order
— Methodorder(t)
Returns the order k
of the knot set t
.
CompactBases.orthogonality
— Functionorthogonality(B)
Returns the BasisOrthogonality
trait of the basis B
.
Examples
julia> orthogonality(BSpline(LinearKnotSet(7, 0, 1, 11)))
CompactBases.NonOrthogonal()
julia> orthogonality(FEDVR(range(0, stop=1, length=11), 7))
CompactBases.Orthonormal()
julia> orthogonality(StaggeredFiniteDifferences(0.1, 0.3, 0.1, 10.0))
CompactBases.Orthonormal()
julia> orthogonality(FiniteDifferences(range(0, stop=1, length=11)))
CompactBases.Orthogonal()
CompactBases.rightmultiplicity
— Methodrightmultiplicity(t)
Return the multiplicity of the right endpoint.
CompactBases.within_interval
— Methodwithin_interval(x, interval)
Return the indices of the elements of x
that lie within the given interval
.
CompactBases.within_support
— Methodwithin_support(x, t, j)
Return the indices of the elements of x
that lie withing the compact support of the j
th basis function (enumerated 1..n
), given the knot set t
. For each index of x
that is covered, the index k
of the interval within which x[i]
falls is also returned.
LinearAlgebra.mul!
— Functionmul!(y, L::DiagonalOperator, x[, α=1, β=0])
Compute the action of the DiagonalOperator
L
on x
and store the result in y
.
CompactBases.@materialize
— Macro@materialize function op(args...)
This macro simplifies the setup of a few functions necessary for the materialization of Applied
objects:
copyto!(dest::DestType, applied_obj::Applied{...,op})
performs the actual materialization ofapplied_obj
into the destination object which has been generated bysimilar
which usually returns a suitable matrixLazyArrays._simplify
which makes use of the above functions
Example
@materialize function *(Ac::MyAdjointBasis,
O::MyOperator,
B::MyBasis)
T -> begin # generates similar
A = parent(Ac)
parent(A) == parent(B) ||
throw(ArgumentError("Incompatible bases"))
# There may be different matrices best representing different
# situations:
if ...
Diagonal(Vector{T}(undef, size(B,1)))
else
Tridiagonal(Vector{T}(undef, size(B,1)-1),
Vector{T}(undef, size(B,1)),
Vector{T}(undef, size(B,1)-1))
end
end
dest::Diagonal{T} -> begin # generate copyto!(dest::Diagonal{T}, ...) where T
dest.diag .= 1
end
dest::Tridiagonal{T} -> begin # generate copyto!(dest::Tridiagonal{T}, ...) where T
dest.dl .= -2
dest.ev .= 1
dest.du .= 3
end
end