Basis Functions
In order to construct a parameterized map, we first need to define a suitable basis.
This manual describes the different one-dimensional basis families currently implemented in TransportMaps.jl. We focus on variants of the (probabilists') Hermite polynomials and recently proposed edge-controlled versions [3], [2].
Probabilistic Hermite Basis
The probabilistic Hermite polynomials form an orthonormal basis with respect to the standard normal density
They satisfy the three-term recurrence
Orthonormal polynomial bases such as the Hermite family are useful for several reasons. Orthogonality with respect to a reference measure makes coefficient estimation stable: projections onto the basis are simple inner products and a truncated series gives the best
We want to visualize the Hermite polynomials, first we load the necessary packages:
using Distributions
using Plots
using TransportMapsThen we construct the basis via HermiteBasis() or HermiteBasis(:none) (:none means no edge control):
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 empirical quantiles of (reference) samples. Alternatively, they can be chosen as the quantiles of the respective reference density. The normalization constant
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")
endSummary
We showcased four basis variants:
Standard (orthonormal) Hermite
Linearized Hermite (piecewise linear tails)
Gaussian-weighted Hermite (exponential damping)
Cubic-spline-weighted Hermite (compact-style smooth damping)
These can be supplied when constructing polynomial transport maps to tune stability, tail behavior, and sparsity characteristics.
This page was generated using Literate.jl.