Basis Functions
In TransportMaps.jl, the basis must match the reference distribution used by the map. The reference should have the same support as the target distribution. If supports do not match, consider an isoprobabilistic transformation before constructing the map.
Supported reference-basis pairings are:
| Reference distribution | Support | Basis family |
|---|---|---|
| Hermite-based bases | ||
| Legendre polynomials | ||
| Shifted Legendre polynomials |
This compatibility is enforced by the map constructor.
Standard Normal Reference
Often, a standard normal reference with the probability density function
is chosen when the reference measure should have an infinite support. In this case, a probabilistic Hermite polynomial basis is a natural choice since they are orthogonal with respect to this density and satisfy the three-term recurrence
with
Orthonormal polynomial bases are useful because projections are stable and truncated expansions give the best
We also support edge-controlled Hermite variants that improve tail behavior while keeping the Gaussian reference [3], [2]:
Probabilistic Hermite basis
Linearized Hermite basis
Gaussian-weighted Hermite basis
Cubic-spline-weighted Hermite basis
Probabilistic Hermite Basis
We first visualize the standard probabilistic Hermite basis.
Construct the basis with HermiteBasis.
basis = HermiteBasis()
z = -3:0.01:3
p1 = plot(xlabel="z", ylabel="Basis function", title="Standard Hermite Basis")
for degree in 0:4
plot!(p1, z, map(x -> basisfunction(basis, degree, x), z), label="degree $degree")
endIf we zoom out, we can see that the tails grow quickly for large
z = -7:0.1:7
p2 = plot(xlabel="z", ylabel="Basis function", title="Standard Hermite Basis")
for degree in 0:4
plot!(p2, z, map(x -> basisfunction(basis, degree, x), z), label="degree $degree")
endLinearized Hermite Basis
Linearized Hermite polynomials (edge-linearized basis) were introduced in [3] to control growth for large
The bounds are chosen here as the 0.01 and 0.99 quantiles of the reference. The normalization follows [3]:
basis = LinearizedHermiteBasis(Normal(), 4, 1)
println("Linearization bounds: ", basis.linearizationbounds)
p3 = plot(xlabel="z", ylabel="Basis function", title="Linearized Hermite Basis")
for degree in 0:4
plot!(p3, z, map(x -> basisfunction(basis, degree, x), z), label="degree $degree")
endLinearization bounds: [-2.326347874040846, 2.326347874040846]Edge-Controlled (Weighted) Hermite Basis: Gaussian Weight
Edge control modifies each Hermite polynomial with a decaying weight to reduce growth in the tails [2]. Using a Gaussian weight gives:
basis = GaussianWeightedHermiteBasis()
p4 = plot(xlabel="z", ylabel="Basis function", title="Gaussian-Weighted Hermite Basis")
for degree in 0:4
plot!(p4, z, map(x -> basisfunction(basis, degree, x), z), label="degree $degree")
endNote
In order to preserve some extrapolation properties, the weights are only applied to polynomials of degree
Edge-Controlled Hermite Basis: Cubic Spline Weight
A cubic spline weight smoothly damps the polynomials outside a radius
basis = CubicSplineHermiteBasis(Normal())
p5 = plot(xlabel="z", ylabel="Basis function", title="Cubic Spline Weighted Hermite Basis")
for degree in 0:4
plot!(p5, z, map(x -> basisfunction(basis, degree, x), z), label="degree $degree")
endUniform[-1, 1] Reference
For a uniform reference on
with
Use this basis with Uniform(-1, 1).
basis = LegendreBasis()
z = -1:0.01:1
p = plot(xlabel="z", ylabel="Basis function", title="Legendre Basis on [-1, 1]")
for degree in 0:4
plot!(p, z, map(x -> basisfunction(basis, degree, x), z), label="degree $degree")
endUniform[0, 1] Reference
For a uniform reference on
Use this basis with Uniform(0, 1).
basis = ShiftedLegendreBasis()
z = 0:0.01:1
p = plot(xlabel="z", ylabel="Basis function", title="Shifted Legendre Basis on [0, 1]")
for degree in 0:4
plot!(p, z, map(x -> basisfunction(basis, degree, x), z), label="degree $degree")
endThis page was generated using Literate.jl.