Bases
Index
TransportMaps.CubicSplineHermiteBasisTransportMaps.GaussianWeightedHermiteBasisTransportMaps.HermiteBasisTransportMaps.LegendreBasisTransportMaps.LinearizedHermiteBasisTransportMaps.MultivariateBasisTransportMaps.ShiftedLegendreBasisTransportMaps.basisfunctionTransportMaps.basisfunction_derivative
Types
TransportMaps.LinearizedHermiteBasis Type
LinearizedHermiteBasisProbabilist Hermite polynomial basis with linearization outside specified bounds.
Fields
linearizationbounds::Vector{Float64}: lower and upper bounds for linearization.normalization::Vector{Float64}: normalization factors for each degree.
Constructors
LinearizedHermiteBasis(; lower=-Inf, upper=Inf, normalization=Float64[]): Explicit constructor with specified bounds and normalization.LinearizedHermiteBasis(max_degree::Int): Construct with default normalization for given maximum degree.LinearizedHermiteBasis(samples::Vector{<:Real}, max_degree::Int, k::Int): Construct bounds from 1st and 99th percentile of samples.LinearizedHermiteBasis(density::Distributions.UnivariateDistribution, max_degree::Int, k::Int): Construct bounds from 1st and 99th percentile of reference density.
TransportMaps.CubicSplineHermiteBasis Type
CubicSplineHermiteBasisProbabilist Hermite polynomial basis with cubic spline edge control.
Fields
radius::Float64: radius of the spline for edge control.
Constructors
CubicSplineHermiteBasis(radius::Float64=3.0): Explicit constructor withradius=3.0CubicSplineHermiteBasis(samples::Vector{<:Real}): Construct radius from from 1st and 99th percentile of samples.CubicSplineHermiteBasis(density::Distributions.UnivariateDistribution): Construct radius from 1st and 99th percentile of reference density.
TransportMaps.GaussianWeightedHermiteBasis Type
GaussianWeightedHermiteBasisProbabilist Hermite polynomial basis with Gaussian weight for edge control.
sourceTransportMaps.LegendreBasis Type
LegendreBasisLegendre polynomial basis, orthogonal with respect to uniform measure on [-1, 1].
sourceTransportMaps.ShiftedLegendreBasis Type
ShiftedLegendreBasisShifted Legendre polynomial basis, orthogonal with respect to uniform measure on [0, 1]. The shifted Legendre polynomials
TransportMaps.MultivariateBasis Type
MultivariateBasis{T<:AbstractPolynomialBasis}A multivariate polynomial basis function constructed as a tensor product of univariate polynomial bases.
Fields
multiindexset::Vector{Int}: Multi-index representing the polynomial degrees for each dimensionunivariatebases::Vector{T}: Vector of univariate basis functions, one for each dimension
Constructors
MultivariateBasis(multiindexset::Vector{Int}, ::Type{T}) where T<:AbstractPolynomialBasis: Create a multivariate basis with the specified multi-index and basis type.MultivariateBasis(multiindexset::Vector{Int}, basistype::AbstractPolynomialBasis): Constructor that accepts a basis instance instead of a type.
Functions
TransportMaps.basisfunction Function
basisfunction(basis::HermiteBasis, αᵢ::Int, zᵢ::Real)Evaluate HermiteBasis with degree αᵢ at zᵢ.
basisfunction(basis::LegendreBasis, αᵢ::Int, zᵢ::Real)Evaluate LegendreBasis with degree αᵢ at zᵢ.
basisfunction(basis::ShiftedLegendreBasis, αᵢ::Int, zᵢ::Real)Evaluate ShiftedLegendreBasis with degree αᵢ at zᵢ ∈ [0,1].
basisfunction(basis::LinearizedHermiteBasis, αᵢ::Int, zᵢ::Real)Evaluate LinearizedHermiteBasis with degree αᵢ at zᵢ.
basisfunction(basis::CubicSplineHermiteBasis, αᵢ::Int, zᵢ::Real)Evaluate CubicSplineHermiteBasis with degree αᵢ at zᵢ.
basisfunction(basis::GaussianWeightedHermiteBasis, αᵢ::Int, zᵢ::Real)Evaluate GaussianWeightedHermiteBasis with degree αᵢ at zᵢ.
TransportMaps.basisfunction_derivative Function
basisfunction_derivative(basis::HermiteBasis, αᵢ::Int, zᵢ::Real)Evaluate derivative of HermiteBasis with degree αᵢ at zᵢ.
basisfunction_derivative(basis::LegendreBasis, αᵢ::Int, zᵢ::Real)Evaluate derivative of LegendreBasis with degree αᵢ at zᵢ.
basisfunction_derivative(basis::ShiftedLegendreBasis, αᵢ::Int, zᵢ::Real)Evaluate derivative of ShiftedLegendreBasis with degree αᵢ at zᵢ ∈ [0,1].
basisfunction_derivative(basis::LinearizedHermiteBasis, αᵢ::Int, zᵢ::Real)Evaluate derivative of LinearizedHermiteBasis with degree αᵢ at zᵢ.
basisfunction_derivative(basis::CubicSplineHermiteBasis, αᵢ::Int, zᵢ::Real)Evaluate derivative of CubicSplineHermiteBasis with degree αᵢ at zᵢ.
basisfunction_derivative(basis::GaussianWeightedHermiteBasis, αᵢ::Int, zᵢ::Real)Evaluate derivative of GaussianWeightedHermiteBasis with degree αᵢ at zᵢ.